CGAL 5.6.1 - Triangulated Surface Mesh Approximation
CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag > Class Template Reference

#include <CGAL/Variational_shape_approximation.h>

Definition

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
TriangleMesha model of FaceListGraph
VertexPointMapa ReadablePropertyMap with boost::graph_traits<TriangleMesh>::vertex_descriptor as key and GeomTraits::Point_3 as value type
ErrorMetricProxya model of ErrorMetricProxy
GeomTraitsa model of Kernel
Concurrency_tagconcurrency tag.
Examples:
Surface_mesh_approximation/vsa_class_interface_example.cpp, and Surface_mesh_approximation/vsa_isotropic_metric_example.cpp.

Types

typedef GeomTraits Geom_traits
 Geometric traits type.
 
typedef ErrorMetricProxy Error_metric
 Error metric for proxy fitting type.
 
typedef Error_metric::Proxy Proxy
 Proxy type.
 
typedef std::array< std::size_t, 3 > Indexed_triangle
 Indexed triangle type.
 

Construction

 Variational_shape_approximation (const TriangleMesh &tm, const VertexPointMap &vpoint_map, const Error_metric &error_metric)
 initializes internal data for the approximation. More...
 

Approximation

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...
 

Refinement Operations

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...
 

Meshing

template<typename NamedParameters >
bool extract_mesh (const NamedParameters &np)
 extracts the output mesh in the form of an indexed triangle set. More...
 

Output

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.
 

Constructor & Destructor Documentation

◆ Variational_shape_approximation()

template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::Variational_shape_approximation ( const TriangleMesh &  tm,
const VertexPointMap &  vpoint_map,
const Error_metric error_metric 
)

initializes internal data for the approximation.

Parameters
tmCGAL TriangleMesh on which approximation operates
vpoint_mapvertex point map of the mesh
error_metrican ErrorMetricProxy object

Member Function Documentation

◆ add_proxies_error_diffusion()

template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
std::size_t CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::add_proxies_error_diffusion ( const std::size_t  nb_proxies)

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_proxiesnumber 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>
std::size_t CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::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.

The re-fitting is performed after each proxy is inserted.

Parameters
nb_proxiesnumber of proxies to be added
nb_iterationsnumber 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>
FT CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::compute_total_error ( )

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 >
bool CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::extract_mesh ( const NamedParameters &  np)

extracts the output mesh in the form of an indexed triangle set.

Template Parameters
NamedParametersa sequence of Named Parameters
Parameters
npan optional sequence of Named Parameters among the ones listed below
Returns
true if the extracted surface mesh is manifold, and false otherwise.
Meshing Named Parameters
  • the chord subdivision ratio threshold to the chord length or average edge length
  • Type: geom_traits::FT
  • Default: 5.0

  • If true, the subdivision_ratio is the ratio of the furthest vertex distance to the chord length, otherwise is the average edge length
  • Type: Boolean
  • Default: false

  • If true, the subdivision_ratio is weighted by dihedral angle
  • Type: Boolean
  • Default: false

  • If true, optimize the anchor locations
  • Type: Boolean
  • Default: true

  • If true, use PCA plane fitting, otherwise use the default area averaged plane parameters
  • Type: Boolean
  • Default: false

◆ find_best_merge()

template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
boost::optional< std::pair<std::size_t, std::size_t> > CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::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.

Note
The best is defined as the minimum merged sum error change (increase or decrease) among all pairs.
Parameters
use_threshold_testif 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 >
std::size_t CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::initialize_seeds ( const NamedParameters &  np)

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
NamedParametersa sequence of Named Parameters
Parameters
npan optional sequence of Named Parameters among the ones listed below
Returns
number of proxies initialized
Seeding Named Parameters
  • the selection of seeding method
  • Type: CGAL::Surface_mesh_approximation::Seeding_method
  • Default: CGAL::Surface_mesh_approximation::HIERARCHICAL

  • the maximum number of proxies used to approximate the input mesh
  • Type: std::size_t
  • Default: num_faces(tm) / 3, used when min_error_drop is also not provided

  • the minimum error drop of the approximation, expressed as the ratio between two iterations of proxy addition
  • Type: geom_traits::FT
  • Default: 0.1, used when max_number_of_proxies is also not provided

  • the number of relaxation iterations interleaved within seeding
  • Type: std::size_t
  • Default: 5

◆ merge()

template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
FT CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::merge ( const std::size_t  px0,
const std::size_t  px1 
)

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
px0the kept proxy index
px1the 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 >
void CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::output ( const NamedParameters &  np) const

outputs approximation results.

Template Parameters
NamedParametersa sequence of Named Parameters
Parameters
npan optional sequence of Named Parameters among the ones listed below
Output Named Parameters
  • a property map to output the proxy index of each face of the input polygon mesh
  • Type: a model of WritablePropertyMap with boost::graph_traits<TriangleMesh>::face_descriptor as key and std::size_t as value type
  • Default: no output operation is performed
  • Extra: A proxy is a set of connected faces which are placed under the same proxy patch (see Figure 76.3)
  • Extra: The proxy-ids are contiguous in range [0, number_of_proxies - 1]

  • an OutputIterator to put proxies in
  • Type: a class model of OutputIterator with CGAL::Surface_mesh_approximation::L21_metric_vector_proxy_no_area_weighting::Proxy value type
  • Default: no output operation is performed

  • an OutputIterator to put anchor points in
  • Type: a class model of OutputIterator with geom_traits::Point_3 value type
  • Default: no output operation is performed

  • an OutputIterator to put indexed triangles in
  • Type: a class model of OutputIterator with std::array<std::size_t, 3> value type
  • Default: no output operation is performed

◆ run()

template<typename TriangleMesh , typename VertexPointMap , typename ErrorMetricProxy = CGAL::Default, typename GeomTraits = CGAL::Default, typename Concurrency_tag = CGAL::Sequential_tag>
FT CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::run ( std::size_t  nb_iterations = 1)

runs the partitioning and fitting processes on the whole surface.

Parameters
nb_iterationsnumber 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>
bool CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::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.

Parameters
cvg_thresholdthe percentage of error change between two successive runs, should be in range (0, 1).
max_iterationsmaximum number of iterations allowed
avg_intervalsize 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>
bool CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::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.

Parameters
px_idxproxy index.
nnumber of split sections.
nb_relaxationsnumber 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>
std::size_t CGAL::Variational_shape_approximation< TriangleMesh, VertexPointMap, ErrorMetricProxy, GeomTraits, Concurrency_tag >::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.

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_proxiesnumber of proxies requested to teleport.
nb_iterationsnumber of re-fitting iterations.
no_threshold_testif 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.