Susan Hert and Stefan Schirra
A subset S ⊆ ℝ3 is convex if for any two points p and q in the set the line segment with endpoints p and q is contained in S. The convex hull of a set S is the smallest convex set containing S. The convex hull of a set of points P ∈ ℝ3 is a convex polytope with vertices in P. A point in P is an extreme point (with respect to P) if it is a vertex of the convex hull of P. A set of points is said to be strongly convex if it consists of only extreme points.
This chapter describes the functions provided in Cgal for producing convex hulls in three dimensions as well as functions for checking if sets of points are strongly convex are not. One can compute the convex hull of a set of points in three dimensions in one of three ways in Cgal: using a static algorithm, using an incremental construction algorithm, or using a triangulation to get a fully dynamic computation.
The function convex_hull_3 provides an implementation of the quickhull algorithm [BDH96] for three dimensions . There are two versions of this function available, one that can be used when it is known that the output will be a polyhedron (i.e., there are more than three points and they are not all collinear) and one that handles all degenerate cases and returns a CGAL::Object, which may be a point, a segment, a triangle, or a polyhedron. Both versions accept a range of input iterators defining the set of points whose convex hull is to be computed and a traits class defining the geometric types and predicates used in computing the hull.
The function convex_hull_3 is parameterized by a traits class, which specifies the types and geometric primitives to be used in the computation. The default for this traits class is Convex_hull_traits_3 .
The function is_strongly_convex_3 implements the algorithm of Mehlhorn et al. [MNS+96] to determine if the vertices of a given polytope constitute a strongly convex point set or not. This function is used in postcondition testing for convex_hull_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/Convex_hull_traits_3.h> #include <CGAL/convex_hull_3.h> #include <vector> typedef CGAL::Exact_predicates_inexact_constructions_kernel 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 if (CGAL::object_cast<Segment_3>(&ch_object) ) std::cout << "convex hull is a segment " << std::endl; else if (CGAL::object_cast<Polyhedron_3>(&ch_object) ) std::cout << "convex hull is a polyhedron " << std::endl; else std::cout << "convex hull error!" << std::endl; return 0; }
The function convex_hull_incremental_3 provides an interface similar to convex_hull_3 for the d-dimensional incremental construction algorithm [CMS93]. implemented by the class CGAL::Convex_hull_d<R> that is specialized to three dimensions. This function accepts an iterator range over a set of input points and returns a polyhedron, but it does not have a traits class in its interface. It uses the kernel class Kernel used in the polyhedron type to define an instance of the adapter traits class CGAL::Convex_hull_d_traits_3<Kernel>.
In most cases, the function convex_hull_3 will be faster than convex_hull_incremental_3. The latter is provided mainly for comparison purposes.
To use the full functionality available with the d-dimensional class CGAL::Convex_hull_d<R> in three dimensions (e.g., the ability to insert new points and to query if a point lies in the convex hull or not), you can instantiate the class CGAL::Convex_hull_d<K> with the adapter traits class CGAL::Convex_hull_d_traits_3<K>, as shown in the following example.
File: demo/Convex_hull_3/incremental_hull_3_demo.cpp
#include <CGAL/Homogeneous.h> #include <CGAL/point_generators_3.h> #include <CGAL/Convex_hull_d.h> #include <CGAL/Convex_hull_d_traits_3.h> #include <CGAL/Convex_hull_d_to_polyhedron_3.h> #include <CGAL/Polyhedron_3.h> #include <CGAL/algorithm.h> #include <CGAL/IO/Geomview_stream.h> #include <CGAL/IO/Polyhedron_geomview_ostream.h> #include <vector> #include <cassert> #ifdef CGAL_USE_GEOMVIEW #ifdef CGAL_USE_LEDA #include <CGAL/leda_integer.h> typedef leda_integer RT; #else #ifdef CGAL_USE_GMP #include <CGAL/Gmpz.h> typedef CGAL::Gmpz RT; #else // NOTE: the choice of double here for a number type may cause problems // for degenerate point sets #include <CGAL/double.h> typedef double RT; #endif #endif typedef CGAL::Homogeneous<RT> K; typedef K::Point_3 Point_3; typedef CGAL::Polyhedron_3< K> Polyhedron_3; typedef CGAL::Convex_hull_d_traits_3<K> Hull_traits_3; typedef CGAL::Convex_hull_d< Hull_traits_3 > Convex_hull_3; typedef CGAL::Creator_uniform_3<double, Point_3> Creator; int main () { Convex_hull_3 CH(3); // create instance of the class with dimension == 3 // generate 250 points randomly on a sphere of radius 100 // and insert them into the convex hull CGAL::Random_points_in_sphere_3<Point_3, Creator> gen(100); for (int i = 0; i < 250 ; i++, ++gen) CH.insert(*gen); assert(CH.is_valid()); // define polyhedron to hold convex hull and create it Polyhedron_3 P; CGAL::convex_hull_d_to_polyhedron_3(CH,P); // display polyhedron in a geomview window CGAL::Geomview_stream geomview; geomview << CGAL::RED; geomview << P; std::cout << "Press any key to end the program: "; std::cout.flush(); char ch; std::cin.get(ch); return 0; } #else int main() { std::cerr << "This demo requires geomview, which is not present on this platform\n"; return 0; } #endif
Fully dynamic maintenance of a convex hull can be achieved by using the class CGAL::Delaunay_triangulation_3. This class supports insertion and removal of points (i.e., vertices of the triangulation) and the convex hull edges are simply the finite edges of infinite faces. The following example illustrates the dynamic construction of a convex hull. First, random points from a sphere of a certain radius are generated and are inserted into a triangulation. Then the number of points of the convex hull are obtained by counting the number of triangulation vertices incident to the infinite vertex. Some of the points are removed and then the number of points remaining on the hull are determined. Notice that the vertices incident to the infinite vertex of the triangulation are on the convex hull but it may be that not all of them are vertices of the hull.
File: examples/Convex_hull_3/dynamic_hull_3.cpp
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h> #include <CGAL/point_generators_3.h> #include <CGAL/Delaunay_triangulation_3.h> #include <CGAL/algorithm.h> #include <list> typedef CGAL::Exact_predicates_inexact_constructions_kernel K; typedef K::Point_3 Point_3; typedef CGAL::Delaunay_triangulation_3<K> Delaunay; typedef Delaunay::Vertex_handle Vertex_handle; int main() { CGAL::Random_points_in_sphere_3<Point_3> gen(100.0); std::list<Point_3> points; // generate 250 points randomly on a sphere of radius 100.0 // and insert them into the triangulation CGAL::copy_n(gen, 250, std::back_inserter(points) ); Delaunay T; T.insert(points.begin(), points.end()); std::list<Vertex_handle> vertices; T.incident_vertices(T.infinite_vertex(), std::back_inserter(vertices)); std::cout << "This convex hull of the 250 points has " << vertices.size() << " points on it." << std::endl; // remove 25 of the input points std::list<Vertex_handle>::iterator v_set_it = vertices.begin(); for (int i = 0; i < 25; i++) { T.remove(*v_set_it); v_set_it++; } vertices.clear(); T.incident_vertices(T.infinite_vertex(), std::back_inserter(vertices)); std::cout << "After removal of 25 points, there are " << vertices.size() << " points on the convex hull." << std::endl; return 0; }