-
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_voxelsrounded down. Cannot be smaller than 512. - Type: unsigned int
- Default: 1,000,000
|
CGAL 6.2 - Convex Decomposition of Polyhedra
|
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. | |
| 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.
| FaceGraph | a model of HalfedgeListGraph, and FaceListGraph |
| OutputIterator | must be an output iterator accepting variables of type std::pair<std::vector<geom_traits::Point_3>, std::vector<std::array<unsigned int, 3> > >. |
| NamedParameters | a sequence of Named Parameters |
| tmesh | the input triangle mesh to approximate by convex volumes |
| out_volumes | output iterator into which convex volumes are recorded |
| np | an optional sequence of Named Parameters among the ones listed below |
| |
| |
| |
| |
| |
| |
| |
| |
|
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.tmesh is a triangle mesh. tmesh is not self-intersecting. tmesh).CGAL::alpha_wrap_3() to generate an enclosing of the input mesh to apply this function onto.CGAL::Polygon_mesh_processing::polygon_soup_to_polygon_mesh() | 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.
| NefPolyhedron_3 | an object of type Nef_polyhedron_3. |
N is bounded. Otherwise, the outer volume is ignored.N is non-convex, it is modified to represent the convex decomposition. If N is convex, it is not modified.CGAL::Nef_polyhedron_3<Traits>