Function

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 not computed.
Precondition: : There are at least four points in the range [first, last) not all of which are collinear.

template <class InputIterator, 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 not 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 . 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

For both versions, if the kernel R of the points determined by InputIterator::value_type 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

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/quickhull_3.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/point_generators_3.h>
#include <CGAL/algorithm.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/convex_hull_3.h>
#include <vector>


typedef CGAL::Exact_predicates_inexact_constructions_kernel  K;
typedef CGAL::Polyhedron_3<K>                     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;

//a functor computing the plane containing a triangular facet
struct Plane_from_facet {
  Polyhedron_3::Plane_3 operator()(Polyhedron_3::Facet& f) {
      Polyhedron_3::Halfedge_handle h = f.halfedge();
      return Polyhedron_3::Plane_3( h->vertex()->point(),
                                    h->next()->vertex()->point(),
                                    h->opposite()->vertex()->point());
  }
};


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::cpp0x::copy_n( gen, 250, std::back_inserter(points) );

  // define polyhedron to hold convex hull
  Polyhedron_3 poly;
  
  // compute convex hull of non-collinear points
  CGAL::convex_hull_3(points.begin(), points.end(), poly);

  std::cout << "The convex hull contains " << poly.size_of_vertices() << " vertices" << std::endl;
  
  // assign a plane equation to each polyhedron facet using functor Plane_from_facet
  std::transform( poly.facets_begin(), poly.facets_end(), poly.planes_begin(),Plane_from_facet());


  return 0;
}