CGAL 5.2.1 - Triangulated Surface Mesh Skeletonization
CGAL::Mean_curvature_flow_skeletonization< TriangleMesh, Traits_, VertexPointMap_, SolverTraits_ > Class Template Reference

#include <CGAL/Mean_curvature_flow_skeletonization.h>

Definition

template<class TriangleMesh, class Traits_ = Default, class VertexPointMap_ = Default, class SolverTraits_ = Default>
class CGAL::Mean_curvature_flow_skeletonization< TriangleMesh, Traits_, VertexPointMap_, SolverTraits_ >

Function object that enables to extract the mean curvature flow skeleton of a triangulated surface mesh.

The algorithm used takes as input a triangulated surface mesh and iteratively contracts the surface mesh following the mean curvature flow [1]. The intermediate contracted surface mesh is called the meso-skeleton. After each iteration, the meso-skeleton is locally remeshed using angle split and edge contraction. The process ends when the modification of the meso-skeleton between two iterations is small.

Template Parameters
TriangleMesha model of FaceListGraph
Traitsa model of MeanCurvatureSkeletonizationTraits
Default:
boost::property_traits<
boost::property_map<TriangleMesh, CGAL::vertex_point_t>::type
>::value_type
VertexPointMapa model of ReadWritePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key and Traits::Point_3 as value type.
Default:
boost::property_map<TriangleMesh, CGAL::vertex_point_t>::const_type.
SolverTraits_a model of NormalEquationSparseLinearAlgebraTraits_d.
Default: If Eigen 3.2 (or greater) is available and CGAL_EIGEN3_ENABLED is defined, then an overload of Eigen_solver_traits is provided as default parameter:
Eigen::SimplicialLDLT< Eigen::SparseMatrix<double> >
>
Examples:
Surface_mesh_skeletonization/MCF_Skeleton_example.cpp, Surface_mesh_skeletonization/MCF_Skeleton_sm_example.cpp, Surface_mesh_skeletonization/segmentation_example.cpp, Surface_mesh_skeletonization/simple_mcfskel_example.cpp, and Surface_mesh_skeletonization/simple_mcfskel_sm_example.cpp.

Types

typedef boost::adjacency_list< boost::vecS, boost::vecS, boost::undirectedS, Vmap > Skeleton
 The graph type representing the skeleton. The vertex property Vmap is a struct with a member point of type Traits::Point_3 and a member vertices of type std::vector<boost::graph_traits<TriangleMesh>::vertex_descriptor>. See the boost documentation page for more details.
 

Constructor

 Mean_curvature_flow_skeletonization (const TriangleMesh &tmesh, VertexPointMap vertex_point_map=get(CGAL::vertex_point, tmesh), Traits traits=Traits())
 The constructor of a skeletonization object. More...
 

Local Remeshing Parameters

double max_triangle_angle () const
 During the local remeshing step, a triangle will be split if it has an angle larger than max_triangle_angle().
 
double min_edge_length () const
 During the local remeshing step, an edge will be collapse if it is length is less than min_edge_length().
 
void set_max_triangle_angle (double value)
 set function for max_triangle_angle()
 
void set_min_edge_length (double value)
 set function for min_edge_length()
 

Algorithm Termination Parameters

std::size_t max_iterations () const
 Maximum number of iterations performed by contract_until_convergence().
 
double area_variation_factor () const
 The convergence is considered to be reached if the variation of the area of the meso-skeleton after one iteration is smaller than area_variation_factor()*original_area where original_area is the area of the input triangle mesh.
 
void set_max_iterations (std::size_t value)
 set function for max_iterations()
 
void set_area_variation_factor (double value)
 set function for area_variation_factor()
 

Vertex Motion Parameters

double quality_speed_tradeoff () const
 This is an advanced function.
Advanced
Controls the velocity of movement and approximation quality: decreasing this value makes the mean curvature flow based contraction converge faster, but results in a skeleton of lower quality. This parameter corresponds to \( w_H \) in the original publication.
.
 
bool is_medially_centered () const
 If true, the meso-skeleton placement will be attracted by an approximation of the medial axis of the mesh during the contraction steps, so will be the result skeleton.
 
double medially_centered_speed_tradeoff () const
 This is an advanced function.
Advanced
Controls the smoothness of the medial approximation: increasing this value results in a (less smooth) skeleton closer to the medial axis, as well as a lower convergence speed. It is only used if is_medially_centered()==true. This parameter corresponds to \( w_M \) in the original publication.
.
 
void set_quality_speed_tradeoff (double value)
 set function for quality_speed_tradeoff()
 
void set_is_medially_centered (bool value)
 set function for is_medially_centered()
 
void set_medially_centered_speed_tradeoff (double value)
 set function for medially_centered_speed_tradeoff()
 

High Level Function

void operator() (Skeleton &skeleton)
 Creates the curve skeleton: the input surface mesh is iteratively contracted until convergence, and then turned into a curve skeleton. More...
 

Low Level Functions

Advanced

The following functions enable the user to run the mean curvature flow skeletonization algorithm step by step.

void contract_geometry ()
 Runs one contraction step following the mean curvature flow.
 
std::size_t collapse_edges ()
 Collapses edges of the meso-skeleton with length less than min_edge_length() and returns the number of edges collapsed.
 
std::size_t split_faces ()
 Splits faces of the meso-skeleton having one angle greater than max_triangle_angle() and returns the number of faces split.
 
std::size_t detect_degeneracies ()
 Prevents degenerate vertices to move during the following contraction steps and returns the number of newly fixed vertices.
 
void contract ()
 Performs subsequent calls to contract_geometry(), collapse_edges(), split_faces() and detect_degeneracies()
 
void contract_until_convergence ()
 Iteratively calls contract() until the change of surface area of the meso-skeleton after one iteration is smaller than area_variation_factor()*original_area where original_area is the area of the input triangle mesh, or if the maximum number of iterations has been reached.
 
void convert_to_skeleton (Skeleton &skeleton)
 Converts the contracted surface mesh to a skeleton curve. More...
 

Access to the Meso-Skeleton

typedef unspecified_type Meso_skeleton
 When using the low level API it is possible to access the intermediate results of the skeletonization process, called meso-skeleton. It is a triangulated surface mesh which is model of FaceListGraph and HalfedgeListGraph.
 
const Meso_skeletonmeso_skeleton () const
 Reference to the collapsed triangulated surface mesh.
 

Constructor & Destructor Documentation

◆ Mean_curvature_flow_skeletonization()

template<class TriangleMesh , class Traits_ = Default, class VertexPointMap_ = Default, class SolverTraits_ = Default>
CGAL::Mean_curvature_flow_skeletonization< TriangleMesh, Traits_, VertexPointMap_, SolverTraits_ >::Mean_curvature_flow_skeletonization ( const TriangleMesh &  tmesh,
VertexPointMap  vertex_point_map = get(CGAL::vertex_point, tmesh),
Traits  traits = Traits() 
)

The constructor of a skeletonization object.

The algorithm parameters are initialized such that:

Precondition
tmesh is a triangulated surface mesh without borders and has exactly one connected component.
Parameters
tmeshinput triangulated surface mesh.
vertex_point_mapproperty map which associates a point to each vertex of the graph.
traitsan instance of the traits class.

Member Function Documentation

◆ convert_to_skeleton()

template<class TriangleMesh , class Traits_ = Default, class VertexPointMap_ = Default, class SolverTraits_ = Default>
void CGAL::Mean_curvature_flow_skeletonization< TriangleMesh, Traits_, VertexPointMap_, SolverTraits_ >::convert_to_skeleton ( Skeleton skeleton)

Converts the contracted surface mesh to a skeleton curve.

Template Parameters
Skeletonan instantiation of boost::adjacency_list as a data structure for the skeleton curve.
Parameters
skeletongraph that will contain the skeleton of tmesh. It should be empty before passed to the function.

◆ operator()()

template<class TriangleMesh , class Traits_ = Default, class VertexPointMap_ = Default, class SolverTraits_ = Default>
void CGAL::Mean_curvature_flow_skeletonization< TriangleMesh, Traits_, VertexPointMap_, SolverTraits_ >::operator() ( Skeleton skeleton)

Creates the curve skeleton: the input surface mesh is iteratively contracted until convergence, and then turned into a curve skeleton.

This is equivalent to calling contract_until_convergence() and convert_to_skeleton().

Parameters
skeletongraph that will contain the skeleton of the input triangulated surface mesh. For each vertex descriptor vd of skeleton, the corresponding point and the set of input vertices that contracted to vd can be retrieved using skeleton[vd].point and skeleton[vd].vertices respectively.