CGAL::convex_hull_3

Definition

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.

#include <CGAL/convex_hull_3.h>

template <class InputIterator, class Polyhedron_3, class Traits>
void
convex_hull_3 ( InputIterator first,
InputIterator last,
Polyhedron_3& P,
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 and the plane equations of each face are computed.
Precondition: : There are at least four points in the range [first, last) not all of which are collinear.

template <class InputIterator, class Polyhedron_3, class Traits>
void
convex_hull_3 ( InputIterator first,
InputIterator last,
Object& ch_object,
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. When the result is a polyhedron, the plane equations of each face are computed.

Requirements

Both functions require the following:
  1. InputIterator::value_type is equivalent to Traits::Point_3.
  2. Traits is a model of the concept ConvexHullTraits_3 . When it is known that the input points are not all coplanar, the types Traits_xy, Traits_yx, and Traits_yz need not be provided. For the purposes of checking the postcondition that the convex hull is valid, Traits should also be a model of the concept IsStronglyConvexTraits_3.

Both functions have an additional requirement for the polyhedron that is to be constructed. For the first version this is that:

and for the second, it is required that

The default traits class for both versions of convex_hull_3 is Convex_hull_traits_3<R>,with the representation R determined by InputIterator::value_type.

See Also

CGAL::convex_hull_incremental_3
CGAL::ch_eddy
CGAL::convex_hull_2

Implementation

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

Example

The following program computes the convex hull of a set of 250 random points chosen from a sphere of radius 100. It then determines if the resulting hull is a segment or a polyhedron. Notice that the traits class is not necessary in the call to convex_hull_3 but is used in the definition of Polyhedron_3.

// file: examples/Convex_hull_3/ch_quickhull_3_ex.C

#include <CGAL/Homogeneous.h>
#include <CGAL/point_generators_3.h>
#include <CGAL/copy_n.h>
#include <CGAL/Convex_hull_traits_3.h>
#include <CGAL/convex_hull_3.h>
#include <vector>

#ifdef CGAL_USE_GMP
#include <CGAL/Gmpz.h>
typedef CGAL::Gmpz RT;
#else
#include <CGAL/MP_Float.h>
typedef CGAL::MP_Float RT;
#endif

typedef CGAL::Homogeneous<RT>                     K;
typedef CGAL::Convex_hull_traits_3<K>             Traits;
typedef Traits::Polyhedron_3                      Polyhedron_3;
typedef K::Segment_3                              Segment_3;

// define point creator 
typedef K::Point_3                                Point_3;
typedef CGAL::Creator_uniform_3<double, Point_3>  PointCreator;


int main()
{
  CGAL::Random_points_in_sphere_3<Point_3, PointCreator> gen(100.0);

  // generate 250 points randomly on a sphere of radius 100.0 
  // and copy them to a vector 
  std::vector<Point_3> points;
  CGAL::copy_n( gen, 250, std::back_inserter(points) );
  
  // define object to hold convex hull 
  CGAL::Object ch_object;

  // compute convex hull 
  CGAL::convex_hull_3(points.begin(), points.end(), ch_object);

  // determine what kind of object it is
  Segment_3 segment;
  Polyhedron_3 polyhedron;
  if ( CGAL::assign(segment, ch_object) )
     std::cout << "convex hull is a segment " << std::endl;
  else if ( CGAL::assign (polyhedron, ch_object) )
     std::cout << "convex hull is a polyhedron " << std::endl;
  else
     std::cout << "convex hull error!" << std::endl;

  return 0;
}