CGAL 4.9 - 3D Convex Hulls
|
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) |
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.
P
.first
, last
) not all of which are collinear.InputIterator | must be an input iterator with a value type equivalent to Traits::Point_3 . |
Polyhedron_3 | must be a model of ConvexHullPolyhedron_3 . |
Traits | must 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.
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>
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.
P
in case the result is a polyhedron.InputIterator | 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 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>
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
.
P
.T.dimension()
==3.Triangulation | is a CGAL 3D triangulation. |
Polyhedron | is an instantiation of CGAL::Polyhedron_3<Traits> . |
convex_hull_3
#include <CGAL/convex_hull_3_to_polyhedron_3.h>
void CGAL::convex_hull_incremental_3 | ( | InputIterator | first, |
InputIterator | beyond, | ||
Polyhedron & | P, | ||
bool | test_correctness = false |
||
) |
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.
Delaunay_triangulation_3
together with convex_hull_3_to_polyhedron_3()
is highly recommended.InputIterator | must be an input iterator where the value type must be Polyhedron::Traits::Point . |
Polyhedron | must provide a type Polyhedron::Traits that defines the following types
|
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>
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.
PlaneIterator
and the point type of origin
must come from the same CGAL Kernel.origin
is inside the intersection of halfspaces defined by the range [begin, end)
. PlaneIterator | must be an input iterator where the value type must be Polyhedron::Traits::Plane_3 |
Polyhedron | must be a model of ConvexHullPolyhedron_3 . |
#include <CGAL/Convex_hull_3/dual/halfspace_intersection_3.h>
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.
PlaneIterator
and the point type of origin
must come from the same CGAL Kernel.origin
is inside the intersection of halfspaces defined by the range [begin, end)
. PlaneIterator | must be an input iterator where the value type must be Polyhedron::Traits::Plane_3 |
Polyhedron | must be a model of ConvexHullPolyhedron_3 . |
Traits | must be a model of the concept ConvexHullTraits_3 . |
halfspace_intersection_3()
#include <CGAL/Convex_hull_3/dual/halfspace_intersection_with_constructions_3.h>