Given two sets , their Minkowski sum, denoted by , is the set . Minkowski sums are used in many applications, such as motion planning and computer-aided design and manufacturing. This package contains functions for computing planar Minkowski sums of two polygons (namely and are two closed polygons in ), and for a polygon and a disc (an operation also known as offsetting or dilating a polygon).
Computing the Minkowski sum of two convex polygons and with and vertices respectively is very easy, as is a convex polygon bounded by copies of the edges, and these edges are sorted by the angle they form with the -axis. As the two input polygons are convex, their edges are already sorted by the angle they form with the -axis. The Minkowski sum can therefore be computed in time, by starting from two bottommost vertices in and in and performing ``merge sort'' on the edges.
If the polygons are not convex, it is possible to use one of the following approaches:
This approach relies on a decomposition strategy that computes the convex decomposition of the input polygons and its performance depends on the quality of the decomposition.
The segments of the convolution form a number of closed (not necessarily simple) polygonal curves called convolution cycles. The Minkowski sum is the set of points having a non-zero winding number with respect to the cycles of .3 See Figure 15.1 for an illustration.
The number of segments in the convolution of two polygons is usually smaller than the number of segments that constitute the boundaries of the sub-sums when using the decomposition approach. As both approaches construct the arrangement of these segments and extract the sum from this arrangement, computing Minkowski sum using the convolution approach usually generates a smaller intermediate arrangement, hence it is faster and consumes less space.
The function minkowski_sum_2 (P, Q) accepts two simple polygons and , represented using the Polygon_2<Kernel,Container> class-template and uses the convolution method in order to compute and return their Minkowski sum .
As the input polygons may not be convex, their Minkowski sum may not be simply connected and contain polygonal holes; see for example Figure 15.1. is therefore an instance of the Polygon_with_holes_2<Kernel,Container> class-template, defined in the Boolean Set-Operations package: The outer boundary of is a polygon that can be accessed using S.outer_boundary(), and its polygonal holes are given by the range [S.holes_begin(), S.holes_end()) (where contains S.number_of_holes() holes in its interior).
The following example program constructs the Minkowski sum of two triangles, as depicted in Figure 15.2. The result in this case is a convex hexagon. This program, as other example programs in this chapter, includes the auxiliary header file ms_rational_nt.h which defines Number_type as either Gmpq or Quotient<MP_Float>, depending on whether the GMP library is installed or not. The file print_util.h includes auxiliary functions for printing polygons.
File: examples/Minkowski_sum_2/sum_triangles.cpp
#include "ms_rational_nt.h" #include <CGAL/Cartesian.h> #include <CGAL/minkowski_sum_2.h> #include <iostream> #include "print_utils.h" typedef CGAL::Cartesian<Number_type> Kernel; typedef Kernel::Point_2 Point_2; typedef CGAL::Polygon_2<Kernel> Polygon_2; typedef CGAL::Polygon_with_holes_2<Kernel> Polygon_with_holes_2; int main () { // Construct the first polygon (a triangle). Polygon_2 P; P.push_back (Point_2 (0, 0)); P.push_back (Point_2 (6, 0)); P.push_back (Point_2 (3, 5)); // Construct the second polygon (a triangle). Polygon_2 Q; Q.push_back (Point_2 (0, 0)); Q.push_back (Point_2 (2, -2)); Q.push_back (Point_2 (2, 2)); // Compute the Minkowski sum. Polygon_with_holes_2 sum = minkowski_sum_2 (P, Q); CGAL_assertion (sum.number_of_holes() == 0); std::cout << "P = "; print_polygon (P); std::cout << "Q = "; print_polygon (Q); std::cout << "P (+) Q = "; print_polygon (sum.outer_boundary()); return (0); }
In the following program we compute the Minkowski sum of two polygons that are read from an input file. In this case, the sum is not simple and contains four holes, as illustrated in Figure 15.3. Note that this example uses the predefined CGAL kernel with exact constructions. In general, instantiating polygons with this kernel yields the fastest running times for Minkowski-sum computations.
File: examples/Minkowski_sum_2/sum_with_holes.cpp
#include <CGAL/Exact_predicates_exact_constructions_kernel.h> #include <CGAL/minkowski_sum_2.h> #include <iostream> #include <fstream> #include "print_utils.h" // instead of //typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel; // workaround for VC++ struct Kernel : public CGAL::Exact_predicates_exact_constructions_kernel {}; typedef Kernel::Point_2 Point_2; typedef CGAL::Polygon_2<Kernel> Polygon_2; typedef CGAL::Polygon_with_holes_2<Kernel> Polygon_with_holes_2; int main () { // Open the input file. std::ifstream in_file ("rooms_star.dat"); if (! in_file.is_open()) { std::cerr << "Failed to open the input file." << std::endl; return (1); } // Read the two polygons from the file and compute their Minkowski sum. Polygon_2 P, Q; in_file >> P >> Q; in_file.close(); // Compute and print the Minkowski sum. Polygon_with_holes_2 sum = minkowski_sum_2 (P, Q); std::cout << "P (+) Q = "; print_polygon_with_holes (sum); return (0); }
In order to compute Minkowski sums using the decomposition method, it is possible to call the function minkowski_sum_2 (P, Q, decomp), where decomp is an instance of a class that models the concept PolygonConvexDecomposition_2, namely it should provide a method named decompose() that receives a planar polygon and returns a range of convex polygons that represents its convex decomposition.
The Minkowski-sum package includes four models of the concept PolygonConvexDecomposition_2. The first three are classes that wrap the decomposition functions included in the Planar Polygon Partitioning package, while the fourth is an implementation of a decomposition algorithm introduced in [AFH02]. The convex decompositions that it creates usually yield efficient running times for Minkowski sum computations:
The following example demonstrates the computation of the Minkowski sum of the same input polygons as used in Minkowski_sum_2/sum_with_holes.cpp (as depicted in Figure 15.3), using the small-side angle-bisector decomposition strategy:
File: examples/Minkowski_sum_2/sum_by_decomposition.cpp
#include <CGAL/Exact_predicates_exact_constructions_kernel.h> #include <CGAL/minkowski_sum_2.h> #include <CGAL/Small_side_angle_bisector_decomposition_2.h> #include <iostream> #include <fstream> #include "print_utils.h" // instead of //typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel; // workaround for VC++ struct Kernel : public CGAL::Exact_predicates_exact_constructions_kernel {}; typedef Kernel::Point_2 Point_2; typedef CGAL::Polygon_2<Kernel> Polygon_2; typedef CGAL::Polygon_with_holes_2<Kernel> Polygon_with_holes_2; int main () { // Open the input file. std::ifstream in_file ("rooms_star.dat"); if (! in_file.is_open()) { std::cerr << "Failed to open the input file." << std::endl; return (1); } // Read the two polygons from the file and compute their Minkowski sum. Polygon_2 P, Q; in_file >> P >> Q; in_file.close(); // Compute the Minkowski sum using the decomposition approach. CGAL::Small_side_angle_bisector_decomposition_2<Kernel> ssab_decomp; Polygon_with_holes_2 sum = minkowski_sum_2 (P, Q, ssab_decomp); std::cout << "P (+) Q = "; print_polygon_with_holes (sum); return (0); }
The operation of computing the Minkowski sum of a polygon with , a disc of radius centered at the origin, is widely known as offsetting the polygon by a radius .
(a) | (b) | (c) |
Let be the polygon vertices, ordered in a counterclockwise orientation around its interior. If is a convex polygon the offset is easily computed by shifting each polygon edge by away from the polygon, namely to the right side of the edge. As a result we obtain a collection of disconnected offset edges. Each pair of adjacent offset edges, induced by and , are connected by a circular arc of radius , whose supporting circle is centered at . The angle that defines such a circular arc equals ; see Figure 15.4(a) for an illustration. The running time of this simple process is of course linear with respect to the size of the polygon.
If is not convex, its offset can be obtained by decomposing it to convex sub-polygons such that , computing the offset of each sub-polygon and finally calculating the union of these sub-offsets (see Figure 15.4(b)). However, as was the case with the Minkowski sum of a pair of polygons, here it is also more efficient to compute the convolution cycle of the polygon with the disc ,4 which can be constructed by applying the process described in the previous paragraph. The only difference is that a circular arc induced by a reflex vertex is defined by an angle ; see Figure 15.4(c) for an illustration. We finally compute the winding numbers of the faces of the arrangement induced by the convolution cycle to obtain the offset polygon.
Let us assume that we are given a counterclockwise-oriented polygon , whose vertices all have rational coordinates (i.e., for each vertex we have , and we wish to compute its Minkowski sum with a disc of radius , where is also a rational number. The boundary of this sum is comprised of line segments and circular arcs, where:
Thus, the linear offset edges are segments of curves of an algebraic curve of degree 2 (a conic curve) with rational coefficients. This curve is actually a pair of the parallel lines . In Section 15.3.2 we use this representation to construct the offset polygon in an exact manner using the traits class for conic arcs.
The class-template Gps_circle_segment_traits_2<Kernel>, included in the Boolean Set-Operations package is used for representing general polygons whose edges are circular arcs or line segments, and for applying set operations (e.g. intersection, union, etc.) on such general polygon. It should be instantiated with a geometric kernel that employs exact rational arithmetic, such that the curves that comprise the polygon edges should be arcs of circles with rational coefficients or segments of lines with rational coefficients. As in our case the line segments do not satisfy this requirement, we apply a simple approximation scheme, such that each irrational line segment is approximated by two rational segments:
The function approximated_offset_2 (P, r, eps) accepts a polygon , an offset radius and . It constructs an approximation for the offset of by the radius using the procedure explained above. Furthermore, it is guaranteed that the approximation error, namely the distance of the point from is bounded by . Using this function, we are capable of working with the Gps_circle_segment_traits_2<Kernel> class, which considerably speeds up the (approximate) construction of the offset polygon and the application of set operations on such polygons. The function returns an object of the nested type Gps_circle_segment_traits_2<Kernel>::Polygon_with_holes_2 representing the approximated offset polygon (recall that if is not convex, its offset may not be simple and may contain holes, whose boundary is also comprised of line segments and circular arcs).
The following example demonstrates the construction of an approximated offset of a non-convex polygon, as depicted in Figure 15.6. Note that we use a geometric kernel parameterized with a filtered rational number-type. Using filtering considerably speeds up the construction of the offset.
File: examples/Minkowski_sum_2/approx_offset.cpp
#include "ms_rational_nt.h" #include <CGAL/Lazy_exact_nt.h> #include <CGAL/Cartesian.h> #include <CGAL/approximated_offset_2.h> #include <CGAL/offset_polygon_2.h> #include <CGAL/Timer.h> #include <iostream> typedef CGAL::Lazy_exact_nt<Number_type> Lazy_exact_nt; // instead of //typedef CGAL::Cartesian<Lazy_exact_nt> Kernel; struct Kernel : public CGAL::Cartesian<Lazy_exact_nt> {}; typedef CGAL::Polygon_2<Kernel> Polygon_2; typedef CGAL::Gps_circle_segment_traits_2<Kernel> Gps_traits_2; typedef Gps_traits_2::Polygon_2 Offset_polygon_2; typedef Gps_traits_2::Polygon_with_holes_2 Offset_polygon_with_holes_2; int main () { // Open the input file. std::ifstream in_file ("spiked.dat"); if (! in_file.is_open()) { std::cerr << "Failed to open the input file." << std::endl; return (1); } // Read the input polygon. Polygon_2 P; in_file >> P; in_file.close(); std::cout << "Read an input polygon with " << P.size() << " vertices." << std::endl; // Approximate the offset polygon. const Number_type radius = 5; const double err_bound = 0.00001; Offset_polygon_with_holes_2 offset; CGAL::Timer timer; timer.start(); offset = approximated_offset_2 (P, radius, err_bound); timer.stop(); std::cout << "The offset polygon has " << offset.outer_boundary().size() << " vertices, " << offset.number_of_holes() << " holes." << std::endl; std::cout << "Offset computation took " << timer.time() << " seconds." << std::endl; return (0); }
As we previously mentioned, it is possible to represent offset polygons in an exact manner, if we treat their edges as arcs of conic curves with rational coefficients. The function offset_polygon_2 (traits, P, r) computes the offset of the polygon by a rational radius in an exact manner. The traits parameter should be an instance of an arrangement-traits class that is capable of handling conic arcs in an exact manner; using the Arr_conic_traits_2 class with the number types provided by the CORE library is the preferred option. The function returns an object of the nested type Gps_traits_2<ArrConicTraits>::Polygons_with_holes_2 (see the documentation of the Boolean Set-Operations package for more details on the traits-class adapter Gps_traits_2), which represented the exact offset polygon.
The following example demonstrates the construction of the offset of the same polygon that serves as an input for the example program Minkowski_sum_2/approx_offset.cpp, presented in the previous subsection (see also Figure 15.6). Note that the resulting polygon is smaller than the one generated by the approximated-offset function (recall that each irrational line segment in this case is approximated by two rational line segments), but the offset computation is considerably slower:
File: examples/Minkowski_sum_2/exact_offset.cpp
#include <CGAL/basic.h> #ifndef CGAL_USE_CORE #include <iostream> int main () { std::cout << "Sorry, this example needs CORE ..." << std::endl; return (0); } #else #include <CGAL/Cartesian.h> #include <CGAL/CORE_algebraic_number_traits.h> #include <CGAL/Arr_conic_traits_2.h> #include <CGAL/offset_polygon_2.h> #include <CGAL/Timer.h> #include <iostream> typedef CGAL::CORE_algebraic_number_traits Nt_traits; typedef Nt_traits::Rational Rational; typedef Nt_traits::Algebraic Algebraic; // instead of //typedef CGAL::Cartesian<Rational> Rat_kernel; //typedef CGAL::Cartesian<Algebraic> Alg_kernel; //typedef CGAL::Arr_conic_traits_2<Rat_kernel, // Alg_kernel, // Nt_traits> Conic_traits_2; // workaround for VC++ struct Rat_kernel : public CGAL::Cartesian<Rational> {}; struct Alg_kernel : public CGAL::Cartesian<Algebraic> {}; struct Conic_traits_2 : public CGAL::Arr_conic_traits_2<Rat_kernel, Alg_kernel, Nt_traits> {}; typedef CGAL::Polygon_2<Rat_kernel> Polygon_2; typedef CGAL::Gps_traits_2<Conic_traits_2> Gps_traits_2; typedef Gps_traits_2::Polygon_2 Offset_polygon_2; typedef Gps_traits_2::Polygon_with_holes_2 Offset_polygon_with_holes_2; int main () { // Open the input file. std::ifstream in_file ("spiked.dat"); if (! in_file.is_open()) { std::cerr << "Failed to open the input file." << std::endl; return (1); } // Read the input polygon. Polygon_2 P; in_file >> P; in_file.close(); std::cout << "Read an input polygon with " << P.size() << " vertices." << std::endl; // Compute the offset polygon. Conic_traits_2 traits; const Rational radius = 5; Offset_polygon_with_holes_2 offset; CGAL::Timer timer; timer.start(); offset = offset_polygon_2 (P, radius, traits); timer.stop(); std::cout << "The offset polygon has " << offset.outer_boundary().size() << " vertices, " << offset.number_of_holes() << " holes." << std::endl; std::cout << "Offset computation took " << timer.time() << " seconds." << std::endl; return (0); } #endif
advanced |
advanced |
1 | Throughout this chapter, we increment or decrement an index of a vertex modulo the size of the polygon. |
2 | We say that a vector lies between two vectors and if we reach strictly before reaching if we move all three vectors to the origin and rotate counterclockwise. Note that this also covers the case where has the same direction as . |
3 | Informally speaking, the winding number of a point with respect to some planar curve is an integer number counting how many times does wind in a counterclockwise direction around . |
4 | As the disc is convex, it is guaranteed that the convolution curve is comprised of a single cycle. |