#include <CGAL/Variational_shape_approximation.h>
template<typename TriangleMesh, typename VertexPointMap, typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
class CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >
Main class for Variational Shape Approximation algorithm.
It is based on [1]. For simple use cases, the function CGAL::Surface_mesh_approximation::approximate_mesh()
might be sufficient.
- Template Parameters
-
- Examples:
- Surface_mesh_approximation/vsa_class_interface_example.cpp, and Surface_mesh_approximation/vsa_isotropic_metric_example.cpp.
|
template<typename NamedParameters > |
std::size_t | initialize_seeds (const NamedParameters &np) |
| initializes the seeds with both maximum number of proxies and minimum error drop stop criteria. More...
|
|
FT | run (std::size_t nb_iterations=1) |
| runs the partitioning and fitting processes on the whole surface. More...
|
|
bool | run_to_convergence (const FT cvg_threshold, const std::size_t max_iterations=100, std::size_t avg_interval=3) |
| calls run while error decrease is greater than cvg_threshold . More...
|
|
FT | compute_total_error () |
| computes fitting error of current partition to the proxies. More...
|
|
std::size_t | add_to_furthest_proxies (const std::size_t nb_proxies, const std::size_t nb_iterations=5) |
| adds proxies to the worst regions one by one. More...
|
|
std::size_t | add_proxies_error_diffusion (const std::size_t nb_proxies) |
| adds proxies by diffusing fitting error into current partition. More...
|
|
|
std::size_t | teleport_proxies (const std::size_t nb_proxies, const std::size_t nb_iterations=5, const bool no_threshold_test=false) |
| teleports the local minimum to the worst region by combining the merging and adding processes. More...
|
|
FT | merge (const std::size_t px0, const std::size_t px1) |
| merges two specified adjacent regions. More...
|
|
boost::optional< std::pair< std::size_t, std::size_t > > | find_best_merge (const bool use_threshold_test) |
| simulates merging and local re-fitting of all pairs of adjacent proxies and finds the best pair to merge. More...
|
|
bool | split (const std::size_t px_idx, const std::size_t n=2, const std::size_t nb_relaxations=10) |
| splits within a specified proxy area via N-section (by default bisection), other regions are not affected. More...
|
|
|
template<typename NamedParameters > |
bool | extract_mesh (const NamedParameters &np) |
| extracts the output mesh in the form of an indexed triangle set. More...
|
|
|
template<typename NamedParameters > |
void | output (const NamedParameters &np) const |
| outputs approximation results. More...
|
|
std::size_t | number_of_proxies () const |
| returns the number of proxies.
|
|
◆ Variational_shape_approximation()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
initializes internal data for the approximation.
- Parameters
-
tm | CGAL TriangleMesh on which approximation operates |
vpoint_map | vertex point map of the mesh |
error_metric | an ErrorMetricProxy object |
◆ add_proxies_error_diffusion()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
adds proxies by diffusing fitting error into current partition.
Each partition is added with the number of proxies in proportion to its fitting error.
- Parameters
-
nb_proxies | number of proxies to be added |
- Returns
- number of proxies successfully added
◆ add_to_furthest_proxies()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
adds proxies to the worst regions one by one.
The re-fitting is performed after each proxy is inserted.
- Parameters
-
nb_proxies | number of proxies to be added |
nb_iterations | number of re-fitting iterations |
- Returns
- number of proxies added
◆ compute_total_error()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
computes fitting error of current partition to the proxies.
- Returns
- total fitting error
◆ extract_mesh()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
template<typename NamedParameters >
extracts the output mesh in the form of an indexed triangle set.
- Template Parameters
-
- Parameters
-
- Returns
true
if the extracted surface mesh is manifold, and false
otherwise.
- Meshing Named Parameters
subdivision_ratio | chord subdivision ratio threshold to the chord length or average edge length. |
relative_to_chord | set true if the subdivision_ratio is the ratio of the furthest vertex distance to the chord length, or to the average edge length otherwise. |
with_dihedral_angle | set true if subdivision_ratio is weighted by dihedral angle. |
optimize_anchor_location | if set to true , optimize the anchor locations. |
pca_plane | set true if use PCA plane fitting, otherwise use the default area averaged plane parameters. |
◆ find_best_merge()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
simulates merging and local re-fitting of all pairs of adjacent proxies and finds the best pair to merge.
- Note
- The best is defined as the minimum merged sum error change (increase or decrease) among all pairs.
- Parameters
-
use_threshold_test | if true and a best pair of proxies is found, it is returned only if the error change after the merge is lower than the half of the maximum proxy error. |
- Returns
- if the best merge pair is found the optional returned contains the proxy indices, and is empty otherwise.
◆ initialize_seeds()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
template<typename NamedParameters >
initializes the seeds with both maximum number of proxies and minimum error drop stop criteria.
The first criterion met stops the seeding. Parameters out of range are ignored.
- Template Parameters
-
- Parameters
-
- Returns
- number of proxies initialized
- Seeding Named Parameters
-
◆ merge()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
merges two specified adjacent regions.
The overall re-fitting is not performed and the proxy index map is maintained.
- Precondition
- two proxies must be adjacent, and 0 <= px0 < px1 < proxies.size()
- Parameters
-
px0 | the kept proxy index |
px1 | the merged and erased proxy index, proxies with greater indeies are decreased by 1 to fill the gap |
- Returns
- change of error
◆ output()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
template<typename NamedParameters >
outputs approximation results.
- Template Parameters
-
- Parameters
-
- Output Named Parameters
face_proxy_map | a WritePropertyMap with boost::graph_traits<TriangleMesh>::face_descriptor as key and std::size_t as value type. A proxy is a set of connected faces which are placed under the same proxy patch (see Figure 71.3). The proxy-ids are contiguous in range [0, number_of_proxies() - 1] . |
proxies | output iterator over proxies. |
anchors | output iterator over anchor points. |
triangles | output iterator over indexed triangles. |
◆ run()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
runs the partitioning and fitting processes on the whole surface.
- Parameters
-
nb_iterations | number of iterations. |
- Returns
- total fitting error
◆ run_to_convergence()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
calls run
while error decrease is greater than cvg_threshold
.
- Parameters
-
cvg_threshold | the percentage of error change between two successive runs, should be in range (0, 1) . |
max_iterations | maximum number of iterations allowed |
avg_interval | size of error average interval to have smoother convergence curve, if 0 is assigned, 1 is used instead. |
- Returns
true
if converged before hitting the maximum iterations.
◆ split()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
splits within a specified proxy area via N-section (by default bisection), other regions are not affected.
- Parameters
-
px_idx | proxy index. |
n | number of split sections. |
nb_relaxations | number of relaxations within the proxy area px_idx after the split |
- Returns
true
if split succeeds, and false
otherwise.
◆ teleport_proxies()
template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
teleports the local minimum to the worst region by combining the merging and adding processes.
The re-fitting is performed after each teleportation. Here if we specify more than one proxy this means we teleport in a naive iterative fashion.
- Parameters
-
nb_proxies | number of proxies requested to teleport. |
nb_iterations | number of re-fitting iterations. |
no_threshold_test | if true , no check on the approximation error before merging a pair of proxies is done. In other words, find_best_merge(!no_threshold_test) is called. |
- Returns
- number of proxies teleported.