\( \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 4.6.3 - 3D Convex Hulls
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
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 Polyhedron >
void CGAL::halfspace_intersection_3 (PlaneIterator begin, PlaneIterator end, Polyhedron &P, boost::optional< Polyhedron::Vertex::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 Polyhedron , class Traits >
void CGAL::halfspace_intersection_with_constructions_3 (PlaneIterator pbegin, PlaneIterator pend, Polyhedron &P, boost::optional< Polyhedron::Vertex::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 Polyhedron_3 , class Traits >
void CGAL::convex_hull_3 (InputIterator first, InputIterator last, Polyhedron_3 &P, 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 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 InputIterator , class Polyhedron >
void CGAL::convex_hull_incremental_3 (InputIterator first, InputIterator beyond, Polyhedron &P, bool test_correctness=false)
 

Function Documentation

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

computes the convex hull of the set of points in the range [first, last).

The polyhedron P is cleared, then the convex hull is stored in P. Note that the convex hull will be triangulated, that is P will contain only triangular facets.

Attention
This function does not compute the plane equations of the faces of P.
Precondition
There are at least four points in the range [first, last) not all of which are collinear.
Template Parameters
InputIteratormust be an input iterator with a value type equivalent to Traits::Point_3.
Polyhedron_3must be a model of ConvexHullPolyhedron_3.
Traitsmust be a model of the concept ConvexHullTraits_3. For the purposes of checking the postcondition that the convex hull is valid, Traits should 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.

See Also
convex_hull_incremental_3()

Implementation

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

#include <CGAL/convex_hull_3.h>

Examples:
Convex_hull_3/quickhull_3.cpp.
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).

The result, which may be a point, a segment, a triangle, or a polyhedron, is stored in ch_object. In the case the result is a polyhedron, the convex hull will be triangulated, that is the polyhedron will contain only triangular facets.

Attention
This function does not compute the plane equations of the faces of P in case the result is a polyhedron.
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 should also be a model of the concept IsStronglyConvexTraits_3. Furthermore, Traits must define a type Polyhedron_3 that is a model of ConvexHullPolyhedron_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.

#include <CGAL/convex_hull_3.h>

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.

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

Attention
This function does not compute the plane equations of the faces of P.
Precondition
T.dimension()==3.
Template Parameters
Triangulationis a CGAL 3D triangulation.
Polyhedronis an instantiation of CGAL::Polyhedron_3<Traits>.
See Also
convex_hull_3

#include <CGAL/convex_hull_3_to_polyhedron_3.h>

Examples:
Convex_hull_3/dynamic_hull_3.cpp.
template<class InputIterator , class Polyhedron >
void CGAL::convex_hull_incremental_3 ( InputIterator  first,
InputIterator  beyond,
Polyhedron &  P,
bool  test_correctness = false 
)
Deprecated:
This function relies on a package deprecated since CGAL 4.6 and thus is also deprecated.

computes the convex hull polyhedron of the three-dimensional points in the range [first,beyond) and assigns it to P. If test_correctness is set to true, the tests described in [3] are used to determine the correctness of the resulting polyhedron.

Attention
This function is provided for completeness and educational purposes. When an efficient incremental implementation is needed, using Delaunay_triangulation_3 together with convex_hull_3_to_polyhedron_3() is highly recommended.
Template Parameters
InputIteratormust be an input iterator where the value type must be Polyhedron::Traits::Point.
Polyhedronmust provide a type Polyhedron::Traits that defines the following types
  • Polyhedron::Traits::R, which is a model of the representation class R required by Convex_hull_d_traits_3<R>
  • Polyhedron::Traits::Point
See Also
convex_hull_3()
Convex_hull_d

Implementation

This function uses the d-dimensional convex hull incremental construction algorithm [2] with d fixed to 3. The algorithm requires \( O(n^2)\) time in the worst case and \( O(n \log n)\) expected time.

#include <CGAL/convex_hull_incremental_3.h>

template<class PlaneIterator , class Polyhedron >
void CGAL::halfspace_intersection_3 ( PlaneIterator  begin,
PlaneIterator  end,
Polyhedron &  P,
boost::optional< Polyhedron::Vertex::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.

The result is stored in the polyhedron P. 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 \) .
The value type of PlaneIterator and the point type of origin must come from the same CGAL Kernel.
Precondition
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 must be Polyhedron::Traits::Plane_3
Polyhedronmust be a model of ConvexHullPolyhedron_3.
See Also
halfspace_intersection_with_constructions_3()

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

Examples:
Convex_hull_3/halfspace_intersection_3.cpp.
template<class PlaneIterator , class Polyhedron , class Traits >
void CGAL::halfspace_intersection_with_constructions_3 ( PlaneIterator  pbegin,
PlaneIterator  pend,
Polyhedron &  P,
boost::optional< Polyhedron::Vertex::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).

The result is stored in the polyhedron P. 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 \) .
The value type of PlaneIterator and the point type of origin must come from the same CGAL Kernel.
Precondition
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 must be Polyhedron::Traits::Plane_3
Polyhedronmust be a model of ConvexHullPolyhedron_3.
Traitsmust be a model of the concept ConvexHullTraits_3.
See Also
halfspace_intersection_3()

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