CGAL 6.2 - Convex Decomposition of Polyhedra
Loading...
Searching...
No Matches
Reference Manual

Peter Hachenberger, Sven Oesau
This packages provides two functions for computing a set of convex volumes that cover a bounded polyhedron. A convex decomposition of Nef Polyhedra into \(O(r^2)\) convex pieces, where \( r\) is the number of edges, whose adjacent facets form an angle of more than 180 degrees with respect to the polyhedron's interior. This bound is worst-case optimal. A second method approximates the input mesh with convex volumes. While these convex volumes cover additional space outside of the polyhedron, the computation is fast for any chosen number of convex volumes.
Introduced in: CGAL 3.5
BibTeX: cgal:h-emspe-26a
License: GPL
Windows Demo: CGAL Lab

Classified Reference Pages

Functions

Functions

template<typename NefPolyhedron_3 >
void CGAL::convex_decomposition_3 (NefPolyhedron_3 &N)
 The function convex_decomposition_3() inserts additional facets into the given Nef polyhedron N, such that each bounded marked volume (the outer volume is unbounded) is subdivided into convex pieces.
 
template<typename FaceGraph , typename OutputIterator , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::approximate_convex_decomposition (const FaceGraph &tmesh, OutputIterator out_volumes, const NamedParameters &np=parameters::default_values())
 approximates a closed input mesh by a set of convex volumes, whose union encloses the input mesh.
 

Function Documentation

◆ approximate_convex_decomposition()

template<typename FaceGraph , typename OutputIterator , typename NamedParameters = parameters::Default_named_parameters>
std::size_t CGAL::approximate_convex_decomposition ( const FaceGraph tmesh,
OutputIterator  out_volumes,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/approximate_convex_decomposition.h>

approximates a closed input mesh by a set of convex volumes, whose union encloses the input mesh.

The bounding box of the input mesh is decomposed into regular voxels. Voxels intersecting the input mesh are then labeled as surface, while the remaining ones are labeled outside or inside according to their position with reference to the domain bounded by the closed input mesh. Next, the convex hull of the mesh is hierarchically split until the volume_error threshold is satisfied or the maximum_depth is reached. Finally, a greedy pairwise merging algorithm combines smaller convex volumes until maximum_number_of_convex_volumes is met.

Template Parameters
FaceGrapha model of HalfedgeListGraph, and FaceListGraph
OutputIteratormust be an output iterator accepting variables of type std::pair<std::vector<geom_traits::Point_3>, std::vector<std::array<unsigned int, 3> > >.
NamedParametersa sequence of Named Parameters
Parameters
tmeshthe input triangle mesh to approximate by convex volumes
out_volumesoutput iterator into which convex volumes are recorded
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • gives an upper bound of the number of voxels. The longest bounding box side will have a length of the cubic root of maximum_number_of_voxels rounded down. Cannot be smaller than 512.
  • Type: unsigned int
  • Default: 1,000,000

  • maximum depth of hierarchical splits
  • Type: unsigned int
  • Default: 10

  • refitting of convex volumes after split phase
  • Type: Boolean
  • Default: true

  • maximum number of convex volumes produced by the method
  • Type: unsigned int
  • Default: 16

  • maximum difference in fraction of volumes of the local convex hull with the sum of voxel volumes. If surpassed, the convex hull will be split if the maximum_depth has not been reached yet.
  • Type: double
  • Default: 0.01

  • If true, the local box of a convex hull is split at the concavity along the longest axis of the bounding box. Otherwise, it is split in the middle of the longest axis, which is faster, but less precise.
  • Type: Boolean
  • Default: true

  • a property map associating points to the vertices of tmesh
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<FaceGraph>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, tmesh)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in FaceGraph.

  • an instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL Kernel deduced from the value type of the vertex-point map, using CGAL::Kernel_traits
  • Extra: The geometric traits class must be compatible with the vertex point type.

Returns
the number of convex hulls. Note that this value may be lower than the maximum_number_of_convex_volumes, for example if the specified volume_error is quickly met. Note that the output volumes can be degenerated, especially, when a large maximum_number_of_convex_volumes is chosen with a limited maximum_depth.
Precondition
tmesh is a triangle mesh.
tmesh is not self-intersecting.
CGAL::is_closed(tmesh).
Note
If the input mesh is not closed, we recommand to first use the function CGAL::alpha_wrap_3() to generate an enclosing of the input mesh to apply this function onto.
See also
CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh()
Examples
Convex_decomposition_3/approximate_convex_decomposition.cpp.

◆ convex_decomposition_3()

template<typename NefPolyhedron_3 >
void CGAL::convex_decomposition_3 ( NefPolyhedron_3 &  N)

#include <CGAL/convex_decomposition_3.h>

The function convex_decomposition_3() inserts additional facets into the given Nef polyhedron N, such that each bounded marked volume (the outer volume is unbounded) is subdivided into convex pieces.

The modified polyhedron represents a decomposition into \(O(r^2)\) convex pieces, where \( r\) is the number of edges that have two adjacent facets that span an angle of more than 180 degrees with respect to the interior of the polyhedron.

The worst-case running time of our implementation is \(O(n^2r^4\sqrt[3]{nr^2}\log{(nr)})\), where \( n\) is the complexity of the polyhedron (the complexity of a Nef polyhedron is the sum of its Vertices, Halfedges and SHalfedges) and \( r\) is the number of reflex edges.

Template Parameters
NefPolyhedron_3an object of type Nef_polyhedron_3.
Precondition
The polyhedron N is bounded. Otherwise, the outer volume is ignored.
Postcondition
If the polyhedron N is non-convex, it is modified to represent the convex decomposition. If N is convex, it is not modified.
See also
CGAL::Nef_polyhedron_3<Traits>
Examples
Convex_decomposition_3/list_of_convex_parts.cpp.