cgal/Surface_mesh_segmentation/doc/Surface_Mesh_Segmentation.txt

104 lines
7.9 KiB
Plaintext

namespace CGAL {
/*!
\mainpage Surface Mesh Segmentation
\authors -
\anchor elephant_sdf_partition
\image latex elephant_sdf_partition.png "Elephant model: SDF values & Segmentation" width=12cm
\image html elephant_sdf_partition.png "Elephant model: SDF values & Segmentation (Going to be changed-improved)"
# Introduction #
Mesh segmentation is the process of partitioning a mesh into smaller and meaningful sub-meshes. Main areas which might utilize mesh segmentation include, but not limited to, modelling, rigging and texturing, shape-retrieval, and deformation.
The algorithm presented here is based on [reference to paper] which introduces Shape Diameter Function that represents a connection between the surface and its volume. Being a pose-invariant and volume-based feature, SDF can be computed for each facet on the surface. Using SDF values, the mesh is first soft clustered by grouping facets that have close SDF values together. These clusters are then refined by graph-cut algorithm which also considers surface-based features like dihedral-angle and concavity.
From the API perspective, the user can access two main outputs of the algorithm: SDF values and segmentation of the mesh. SDF values might be useful for some other purposes than segmentation such as skeletonization, and part-retrieval.
# Overview of the Segmentation Process #
The segmentation algorithm consists of three major parts: Shape Diameter Function (SDF), soft clustering, and graph-cut for hard clustering.
## Shape Diameter Function ##
First part is computation of Shape Diameter Function (SDF), which represents a connection between the surface and its volume. More precisely, SDF is a scalar-valued function defined on the surface which measures the corresponding local volume.
The main handiness of SDF is being able to distinguish thick and thin parts of the mesh by bringing in a volume-based feature to the surface. Another key feature of SDF is its pose-invariant nature, which means that SDF values for a mesh remains largely uneffected after changes of pose.
Basically, SDF for a surface is computed by processing each facets, which comprise the surface, one by one. For a given facet, SDF value computation begins with casting several rays sampled from a cone which is constructed using centroid of the facet as apex and inward-normal of the facet as axis. Then using casted rays (which intuitively correspond to sampling for local volume), resulting SDF value is calculated by first applying outlier removal and then taking weighted average of ray distances.
After calculating SDF values for each facet, bilateral smoothing (an edge-preserving filtering technique) is applied on SDF values. The purpose of edge-preserving smoothing is keeping SDF edges (fast changes on SDF values) in-place without blurring, since they are good candidates for segment boundaries.
## Soft Clustering ##
Second part consists of applying soft clustering on SDF values. For this purpose, SDF values are first clustered using k-means clustering algorithm which is initialized with k-means++ (an algorithm for choosing random seeds for clusters) and run multiple times. Among these runs, we choose clustering result that has minimum with-in cluster error and use it to initialize expectation maximization algorithm for fitting Gaussian mixture models.
Note that there is no direct relation between number of clusters (parameter for soft clustering) and number of segments (disconnected components). Intuitively, number of clusters represents level of segmentation by clustering facets, which have close SDF values, together without considering their connectivity. For example, for a human model, 2 clusters might result in clustering 4 legs and head of the model into one cluster, and remaining parts into another cluster. However, large number of clusters likely to result in detailed segmentation of the mesh with large number of segments.
At the end of this part, we have soft clustering result which is a matrix that contains probability values for each facet to each cluster. These probability values are used in graph-cut part.
## Graph-Cut ##
For hard clustering, which gives us the final partitioning of the mesh, an energy function that combines probability matrix and surface features is minimized. The aim of energy function is while trying to assign facets to their most probable clusters, it is also trying to minimize cost of boundary edges by considering their concavity and dihedral-angle (concave edges with large dihedral angles are considered to be less costly, which means that they are good candidates for segment boundaries).
Energy function that is minimized using alpha-expansion graph cut algorithm is:
<table border="0">
<tr>
<td>
\f$ E(\bar{x}) = \sum\limits_{f \in F} e_1(f, x_f) + \lambda \sum\limits_{ \{f,g\} \in N} e_2(x_f, x_g) \f$
\f$ e_1(f, x_f) = -log(max(P(f|x_f), \epsilon)) \f$
\f$ e_2(x_f, x_g) =
\left \{
\begin{array}{rl}
-log(\theta(f,g)/\pi) &\mbox{ $x_f \ne x_g$} \\
0 &\mbox{ $x_f = x_g$}
\end{array}
\right \} \f$
</td>
<td>
where:
- \f$F\f$ represents set of facets,
- \f$N\f$ represents set of pairs of neighbor facets,
- \f$x_f\f$ is assigned cluster to facet \f$f\f$,
- \f$P(f|x_p)\f$ is probability of assigning facet \f$f\f$ to cluster \f$x_p\f$,
- \f$\theta(f,g)\f$ is dihedral angle between facets \f$f\f$, and \f$g\f$.
</td>
</tr>
</table>
Having facets assigned to corresponding clusters, segments are determined by giving a unique id to each set of connected facets which are placed under same cluster.
# API #
The API provides three free template functions listed below:
- CGAL::sdf_values_computation : computes SDF values.
- CGAL::surface_mesh_segmentation_from_sdf_values : computes mesh segmentation from SDF values.
- CGAL::surface_mesh_segmentation : computes SDF values and mesh segmentation in one go (i.e. combines two functions listed above).
##SDF Computation##
The function CGAL::sdf_values_computation provides an implementation to SDF computation for a given `CGAL Polyhedron`. After computation, following postprocesses 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].
###Examples####
\include Surface_mesh_segmentation/sdf_values_computation_example.cpp
##Surface Mesh Segmentation##
The function CGAL::surface_mesh_segmentation_from_sdf_values implements surface mesh segmentation algorithm. By taking SDF values and segmentation parameters as input, a map which contains a segment-id (between [0, number of segments -1]) for each facet is provided as output. Formally, a segment is set of connected facets which are placed under same cluster after graph-cut.
###Examples###
\include Surface_mesh_segmentation/surface_mesh_segmentation_from_sdf_values_example.cpp
The function CGAL::surface_mesh_segmentation combines two functions explained above 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 same SDF values.
###Examples###
\include Surface_mesh_segmentation/surface_mesh_segmentation_example.cpp
# Extra: (ideas to append) #
Talk about effect of "number of rays" on SDF values, by providing RSME and average error ?
Talk about effect of "number of levels" on segmentation results ? giving a few screen shot.
Talk about pose-invariant feature of segmentation by giving a few screen shots from "benchmark" models (same models with different poses) ?
Adding results of "benchmark" at the end ? Screenshots ?
*/
} /* namespace CGAL */