The function convex_hull_3() computes the convex hull of a given set of three-dimensional points.
Two versions of this function are available. The first can be used when it is known that the result will be a polyhedron and the second when a degenerate hull may also be possible.
Functions
template<class PlaneIterator , class PolygonMesh >
computes robustly the intersection of the halfspaces defined by the planes contained in the range [begin, end) without constructing the dual points. More...
template<class PlaneIterator , class PolygonMesh , class Traits >
computes the convex hull of the set of points in the range [first, last).
The polygon mesh pm is cleared, then the convex hull is stored in pm. Note that the convex hull will be triangulated, that is pm will contain only triangular facets. if the convex hull is a point or a segment, endpoints will be added in pm as isolated vertices.
must be a model of the concept ConvexHullTraits_3. For the purposes of checking the postcondition that the convex hull is valid, Traits must also be a model of the concept IsStronglyConvexTraits_3.
If the kernel R of the points determined by the value type of InputIterator is a kernel with exact predicates but inexact constructions (in practice we check R::Has_filtered_predicates_tag is Tag_true and R::FT is a floating point type), then the default traits class of convex_hull_3() is Convex_hull_traits_3<R>, and R otherwise.
Attention
The user must include the header file of the PolygonMesh type.
Implementation
The algorithm implemented by these functions is the quickhull algorithm of Barnard et al.[1].
computes the convex hull of the points associated to the vertices of g.
The polygon mesh pm is cleared, then the convex hull is stored in pm. Note that the convex hull will be triangulated, that is pm will contain only triangular faces. if the convex hull is a point or a segment, endpoints will be added in pm as isolated vertices.
computes the convex hull of the set of points in the range [first, last).
The result, which may be a point, a segment, a triangle, or a polygon mesh, is stored in ch_object. In the case the result is a polygon mesh, the convex hull will be triangulated, that is the polygon mesh will contain only triangular facets.
must be an input iterator with a value type equivalent to Traits::Point_3.
Traits
must be model of the concept ConvexHullTraits_3. For the purposes of checking the postcondition that the convex hull is valid, Traits must also be a model of the concept IsStronglyConvexTraits_3. Furthermore, Traits must define a type PolygonMesh that is a model of MutableFaceGraph.
If the kernel R of the points determined by the value type of InputIterator is a kernel with exact predicates but inexact constructions (in practice we check R::Has_filtered_predicates_tag is Tag_true and R::FT is a floating point type), then the default traits class of convex_hull_3() is Convex_hull_traits_3<R>, and R otherwise.
Attention
The user must include the header file of the PolygonMesh type.
If the kernel R of the points from InputRange is a kernel with exact predicates but inexact constructions (in practice we check R::Has_filtered_predicates_tag is Tag_true and R::FT is a floating point type), then the default traits class used is Convex_hull_traits_3<R>, and R otherwise.
Parameters
range
the range of input points.
out
an output iterator where the extreme points will be put.
computes robustly the intersection of the halfspaces defined by the planes contained in the range [begin, end) without constructing the dual points.
The result is stored in the polyhedron pm. If origin is given then it must be a point strictly inside the polyhedron. If an interior point is not given then it is computed using the function halfspace_intersection_interior_point_3() based on solving a linear program and thus is slower.
This version does not construct the dual points explicitely but uses a special traits class for the function CGAL::convex_hull_3() to handle predicates on dual points without constructing them.
Halfspaces are considered as lower halfspaces, that is if the plane equation is \( a\, x +b\, y +c\, z + d = 0 \) then the corresponding halfspace is defined by \( a\, x +b\, y +c\, z + d \le 0 \) .
Precondition
The point type of origin and the point type of the vertices of PolygonMesh must come from the same CGAL Kernel.
if provided, origin is inside the intersection of halfspaces defined by the range [begin, end).
The computed intersection must be a bounded convex polyhedron.
Template Parameters
PlaneIterator
must be an input iterator where the value type is a model of the concept Kernel::Plane_3 and this plane type must come from the same kernel as the point type.
computes a point belonging to the intersection of the halfspaces defined by the planes contained in the range [begin, end).
If the intersection is empty, boost::none is returned.
Attention
Halfspaces are considered as lower halfspaces that is to say if the plane's equation is \( a\, x +b\, y +c\, z + d = 0 \) then the corresponding halfspace is defined by \( a\, x +b\, y +c\, z + d \le 0 \) .
Template Parameters
PlaneIterator
must be an input iterator with the value type being a Plane_3 object from CGAL Kernel
computes the intersection of the halfspaces defined by the planes contained in the range [begin, end).
The result is stored in the polyhedron pm. If origin is given then it must be a point strictly inside the polyhedron. If an interior point is not given then it is computed using the function halfspace_intersection_interior_point_3() based on solving a linear program and thus is slower. This version constructs explicitly the dual points using the convex hull algorithm parametrized with the given traits class.
Halfspaces are considered as lower halfspaces, that is if the plane equation is \( a\, x +b\, y +c\, z + d = 0 \) then the corresponding halfspace is defined by \( a\, x +b\, y +c\, z + d \le 0 \) .
Precondition
The value type of PlaneIterator and the point type of origin must come from the same CGAL Kernel.
if provided, origin is inside the intersection of halfspaces defined by the range [begin, end).
The computed intersection must be a bounded convex polyhedron.
Template Parameters
PlaneIterator
must be an input iterator where the value type is a model of the concept Kernel::Plane_3 and this plane type must come from the same kernel as the point type.