CGAL 6.1 - 3D Convex Hulls
Loading...
Searching...
No Matches
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 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).
 
template<class VertexListGraph , class PolygonMesh , class NamedParameters = parameters::Default_named_parameters>
void CGAL::convex_hull_3 (const VertexListGraph &g, PolygonMesh &pm, const NamedParameters &np=parameters::default_values())
 computes the convex hull of the points associated to the vertices of g.
 
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).
 
template<class InputIterator , class PointRange , class PolygonRange class Traits>
void CGAL::convex_hull_3 (InputIterator first, InputIterator last, PointRange &vertices, PolygonRange &faces, const Traits &ch_traits=Default_traits)
 computes the convex hull of the set of points in the range [first, last).
 
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.
 
template<class PlaneIterator , class PolygonMesh >
void CGAL::halfspace_intersection_3 (PlaneIterator begin, PlaneIterator end, PolygonMesh &pm, std::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 > > origin=std::nullopt)
 computes robustly the intersection of the halfspaces defined by the planes contained in the range [begin, end) without constructing the dual points.
 
template<class PlaneIterator , class PolygonMesh , class Traits >
void CGAL::halfspace_intersection_with_constructions_3 (PlaneIterator pbegin, PlaneIterator pend, PolygonMesh &pm, std::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 > > origin=std::nullopt, const Traits &ch_traits=Default_traits)
 computes the intersection of the halfspaces defined by the planes contained in the range [begin, end).
 
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.
 
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).
 
template<class PlaneIterator >
std::optional< typename Kernel_traits< typename std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3CGAL::halfspace_intersection_interior_point_3 (PlaneIterator begin, PlaneIterator end)
 computes a point belonging to the intersection of the halfspaces defined by the planes contained in the range [begin, end).
 

Function Documentation

◆ convex_hull_3() [1/4]

template<class VertexListGraph , class PolygonMesh , class NamedParameters = parameters::Default_named_parameters>
void CGAL::convex_hull_3 ( const VertexListGraph g,
PolygonMesh &  pm,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/convex_hull_3.h>

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.

Template Parameters
VertexListGrapha model of VertexListGraph.
PolygonMeshmust be a model of MutableFaceGraph. an internal property map for CGAL::vertex_point_t must be available.
NamedParametersa sequence of named parameters
Parameters
gthe graph
pmthe PolygonMesh that will contain the convex hull
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of g
  • Type: a class model of ReadablePropertyMap with boost::graph_traits<VertexListGraph>::vertex_descriptor as key type and Point_3 as value type
  • Default: boost::get(CGAL::vertex_point, g)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in VertexListGraph.
Attention
The user must include the header file of the PolygonMesh and VertexListGraph types.

◆ convex_hull_3() [2/4]

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

◆ convex_hull_3() [3/4]

template<class InputIterator , class PointRange , class PolygonRange class Traits>
void CGAL::convex_hull_3 ( InputIterator  first,
InputIterator  last,
PointRange &  vertices,
PolygonRange &  faces,
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 triangle mesh, is stored as an indexed triangle soup in a vector of points and a vector of index triples.

Template Parameters
InputIteratormust be an input iterator with a value type equivalent to Traits::Point_3.
PointRangea model of the concept RandomAccessContainer whose value type is the point type.
PolygonRangea model of the concepts RandomAccessContainer and BackInsertionSequence whose value_type is itself a model of the concepts RandomAccessContainer and BackInsertionSequence whose value_type is an unsigned integer type convertible to std::size_t
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.

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.

◆ convex_hull_3() [4/4]

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

Implementation

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

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

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

◆ 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,
std::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 ,
origin  = std::nullopt 
)

#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 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 explicitly 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
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_interior_point_3()

template<class PlaneIterator >
std::optional< typename Kernel_traits< typename std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 > CGAL::halfspace_intersection_interior_point_3 ( PlaneIterator  begin,
PlaneIterator  end 
)

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

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, std::nullopt 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
PlaneIteratormust be an input iterator with the value type being a Plane_3 object from CGAL Kernel

◆ 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,
std::optional< Kernel_traits< std::iterator_traits< PlaneIterator >::value_type >::Kernel::Point_3 ,
origin  = std::nullopt,
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 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
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()