\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 5.0 - 3D Convex Hulls
Convex Hull Functions

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 >
void CGAL::halfspace_intersection_3 (PlaneIterator begin, PlaneIterator end, PolygonMesh &pm, boost::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 > > origin=boost::none)
 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 >
void CGAL::halfspace_intersection_with_constructions_3 (PlaneIterator pbegin, PlaneIterator pend, PolygonMesh &pm, boost::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 > > origin=boost::none, const Traits &ch_traits=Default_traits)
 computes the intersection of the halfspaces defined by the planes contained in the range [begin, end). More...
 
template<class InputIterator , class PolygonMesh , class Traits >
void CGAL::convex_hull_3 (InputIterator first, InputIterator last, PolygonMesh &pm, const Traits &ch_traits=Default_traits)
 computes the convex hull of the set of points in the range [first, last). More...
 
template<class InputIterator , class Traits >
void CGAL::convex_hull_3 (InputIterator first, InputIterator last, Object &ch_object, const Traits &ch_traits=Default_traits)
 computes the convex hull of the set of points in the range [first, last). More...
 
template<class InputRange , class OutputIterator , class Traits >
OutputIterator CGAL::extreme_points_3 (InputRange range, OutputIterator out, const Traits &traits)
 copies in out the points on the convex hull of the points in range. More...
 
template<class Triangulation , class PolygonMesh >
void CGAL::convex_hull_3_to_face_graph (const Triangulation &T, PolygonMesh &pm)
 fills a polyhedron with the convex hull of a set of 3D points contained in a 3D triangulation of CGAL. More...
 
template<class Triangulation , class Polyhedron >
void CGAL::convex_hull_3_to_polyhedron_3 (const Triangulation &T, Polyhedron &P)
 fills a polyhedron with the convex hull of a set of 3D points contained in a 3D triangulation of CGAL. More...
 
template<class PointPropertyMap , class Base_traits >
Extreme_points_traits_adapter_3< PointPropertyMap, Base_traits > CGAL::make_extreme_points_traits_adapter (const PointPropertyMap &pmap, Base_traits traits)
 Returns Extreme_points_traits_adapter_3<PointPropertyMap, Base_traits>(pmap, traits).
 

Function Documentation

◆ convex_hull_3() [1/2]

template<class InputIterator , class PolygonMesh , class Traits >
void CGAL::convex_hull_3 ( InputIterator  first,
InputIterator  last,
PolygonMesh &  pm,
const Traits &  ch_traits = Default_traits 
)

#include <CGAL/convex_hull_3.h>

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.

Template Parameters
InputIteratormust be an input iterator with a value type equivalent to Traits::Point_3.
PolygonMeshmust be a model of MutableFaceGraph.
Traitsmust 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 Polygon_mesh type.

Implementation

The algorithm implemented by these functions is the quickhull algorithm of Barnard et al. [1].

Examples:
Convex_hull_3/quickhull_3.cpp, Convex_hull_3/quickhull_any_dim_3.cpp, and Convex_hull_3/quickhull_OM_3.cpp.

◆ convex_hull_3() [2/2]

template<class InputIterator , class Traits >
void CGAL::convex_hull_3 ( InputIterator  first,
InputIterator  last,
Object ch_object,
const Traits &  ch_traits = Default_traits 
)

#include <CGAL/convex_hull_3.h>

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.

Template Parameters
InputIteratormust be an input iterator with a value type equivalent to Traits::Point_3.
Traitsmust 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 Polygon_mesh 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 Polygon_mesh type.

◆ convex_hull_3_to_face_graph()

template<class Triangulation , class PolygonMesh >
void CGAL::convex_hull_3_to_face_graph ( const Triangulation &  T,
PolygonMesh &  pm 
)

#include <CGAL/convex_hull_3_to_face_graph.h>

fills a polyhedron with the convex hull of a set of 3D points contained in a 3D triangulation of CGAL.

The polyhedron pm is cleared and the convex hull of the set of 3D points is stored in pm.

Precondition
T.dimension()==3.
Template Parameters
Triangulationmust be a CGAL 3D triangulation
PolygonMeshmust be a model of the concept MutableFaceGraph
See also
convex_hull_3()
link_to_face_graph()
Examples:
Convex_hull_3/dynamic_hull_3.cpp, and Convex_hull_3/dynamic_hull_OM_3.cpp.

◆ convex_hull_3_to_polyhedron_3()

template<class Triangulation , class Polyhedron >
void CGAL::convex_hull_3_to_polyhedron_3 ( const Triangulation &  T,
Polyhedron &  P 
)

#include <CGAL/convex_hull_3_to_polyhedron_3.h>

fills a polyhedron with the convex hull of a set of 3D points contained in a 3D triangulation of CGAL.

The polyhedron P is cleared and the convex hull of the set of 3D points is stored in P.

Deprecated:
since CGAL 4.10. Use convex_hull_3_to_face_graph() instead.
Attention
This function does not compute the plane equations of the faces of P.
This function works only for CGAL::Polyhedron_3<Traits>, and users who want to generate a Surface_mesh or any other model of a FaceGraph may use convex_hull_3_to_face_graph() instead.
Precondition
T.dimension()==3.
Template Parameters
Triangulationis a CGAL 3D triangulation.
Polyhedronis an instantiation of CGAL::Polyhedron_3<Traits>.
See also
convex_hull_3()
link_to_face_graph()

◆ extreme_points_3()

template<class InputRange , class OutputIterator , class Traits >
OutputIterator CGAL::extreme_points_3 ( InputRange  range,
OutputIterator  out,
const Traits &  traits 
)

#include <CGAL/convex_hull_3.h>

copies in out the points on the convex hull of the points in range.

Template Parameters
InputRangea range of Traits::Point_3, model of ConstRange. Its iterator type is InputIterator.
OutputIteratormust be an output iterator where points of type Traits::Point_3 can be put.
Traitsmust be model of the concept ConvexHullTraits_3.

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
rangethe range of input points.
outan output iterator where the extreme points will be put.
traitsan instance of Traits.
See also
CGAL::Extreme_points_traits_adapter_3
Examples:
Convex_hull_3/extreme_indices_3.cpp, and Convex_hull_3/extreme_points_3_sm.cpp.

◆ halfspace_intersection_3()

template<class PlaneIterator , class PolygonMesh >
void CGAL::halfspace_intersection_3 ( PlaneIterator  begin,
PlaneIterator  end,
PolygonMesh &  pm,
boost::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 ,
origin  = boost::none 
)

#include <CGAL/Convex_hull_3/dual/halfspace_intersection_3.h>

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

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 \) .
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
PlaneIteratormust 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.
PolygonMeshmust be a model of MutableFaceGraph.
See also
halfspace_intersection_with_constructions_3()
Examples:
Convex_hull_3/halfspace_intersection_3.cpp.

◆ halfspace_intersection_with_constructions_3()

template<class PlaneIterator , class PolygonMesh , class Traits >
void CGAL::halfspace_intersection_with_constructions_3 ( PlaneIterator  pbegin,
PlaneIterator  pend,
PolygonMesh &  pm,
boost::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 ,
origin  = boost::none,
const Traits &  ch_traits = Default_traits 
)

#include <CGAL/Convex_hull_3/dual/halfspace_intersection_with_constructions_3.h>

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

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 \) .
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
PlaneIteratormust 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.
PolygonMeshmust be a model of MutableFaceGraph.
Traitsmust be a model of the concept ConvexHullTraits_3.
See also
halfspace_intersection_3()