Chapter 22
2D Minkowski Sums

Ron Wein

Table of Contents

22.1 Introduction
22.2 Computing the Minkowski Sum of Two Polygons
   22.2.1   Computing Minkowski Sum using Convolutions
   22.2.2   Decomposition Strategies
22.3 Offsetting a Polygon
   22.3.1   Approximating the Offset with a Guaranteed Error-Bound
   22.3.2   Computing the Exact Offset
   22.3.3   Computing Inner Offsets

22.1   Introduction

Given two sets A,B d, their Minkowski sum, denoted by A B, is the set { a + b  |  a A, b B }. 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 A and B are two closed polygons in 2), and for a polygon and a disc (an operation also known as offsetting or dilating a polygon).

22.2   Computing the Minkowski Sum of Two Polygons

Computing the Minkowski sum of two convex polygons P and Q with m and n vertices respectively is very easy, as P Q is a convex polygon bounded by copies of the m + n edges, and these edges are sorted by the angle they form with the x-axis. As the two input polygons are convex, their edges are already sorted by the angle they form with the x-axis. The Minkowski sum can therefore be computed in O(m + n) time, by starting from two bottommost vertices in P and in Q and performing ``merge sort'' on the edges.

Convolution cycle Convolution cycle
Figure 22.1:  Computing the convolution of a convex polygon and a non-convex polygon (left). The convolution consists of a single self-intersecting cycle, drawn as a sequence of arrows (right). The winding number associated with each face of the arrangement induced by the segments forming the cycle appears in dashed circles. The Minkowski sum of the two polygons is shaded.

If the polygons are not convex, it is possible to use one of the following approaches:

Decomposition:
We decompose P and Q into convex sub-polygons, namely we obtain two sets of convex polygons P1, , Pk and Q1, , Q such that i = 1kPi = P and i = jQj = Q. We then calculate the pairwise sums Sij = Pi Qj using the simple procedure described above, and compute the union P Q = ijSij.

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.

Convolution:
Let us denote the vertices of the input polygons by P = ( p0, , pm-1 ) and Q = ( q0, , qn-1 ). We assume that both P and Q have positive orientations (i.e. their boundaries wind in a counterclockwise order around their interiors) and compute the convolution of the two polygon boundaries. The convolution of these two polygons [GRS83], denoted P * Q, is a collection of line segments of the form [pi + qj, pi+1 + qj],1 where the vector pi pi+1 lies between qj-1 qj and qj qj+1,2 and - symmetrically - of segments of the form [pi + qj, pi + qj+1], where the vector qj qj+1 lies between pi-1 pi and pi pi+1.

The segments of the convolution form a number of closed (not necessarily simple) polygonal curves called convolution cycles. The Minkowski sum P Q is the set of points having a non-zero winding number with respect to the cycles of P * Q.3 See Figure 22.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 Sij 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.

22.2.1   Computing Minkowski Sum using Convolutions

The function minkowski_sum_2 (P, Q) accepts two simple polygons P and Q, represented using the Polygon_2<Kernel,Container> class-template and uses the convolution method in order to compute and return their Minkowski sum S = P Q.

As the input polygons may not be convex, their Minkowski sum may not be simply connected and contain polygonal holes; see for example Figure 22.1. S 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 S 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 S contains S.number_of_holes() holes in its interior).

Minkowski sum of two triangles
Figure 22.2:  Computing the Minkowski sum of two triangles, as done in the example program Minkowski_sum_2/sum_triangles.cpp.

The following example program constructs the Minkowski sum of two triangles, as depicted in Figure 22.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);
}

Minkowski sum of two non-convex polygon that contains holes
Figure 22.3:  Computing the Minkowski sum of two non-convex polygons P and Q, as done in the example programs Minkowski_sum_2/sum_with_holes.cpp and Minkowski_sum_2/sum_by_decomposition.cpp.

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 22.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"

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);
}

22.2.2   Decomposition Strategies

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 22.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"

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);
}

22.3   Offsetting a Polygon

The operation of computing the Minkowski sum P Br of a polygon P with br, a disc of radius r centered at the origin, is widely known as offsetting the polygon P by a radius r.

Computing offsets of polygons Computing offsets of polygons Computing offsets of polygons
(a)(b)(c)
Figure 22.4:  (a) Offsetting a convex polygon. (b) Computing the offset of a non-convex polygon by decomposing it to convex sub-polygons. (c) Offsetting a non-convex polygon by computing its convolution with a disc. The convolution cycle induces an arrangement with three faces, whose winding numbers are shown enclosed in dashed circles.

Let P = ( p0, , pn - 1 ) be the polygon vertices, ordered in a counterclockwise orientation around its interior. If P is a convex polygon the offset is easily computed by shifting each polygon edge by r away from the polygon, namely to the right side of the edge. As a result we obtain a collection of n disconnected offset edges. Each pair of adjacent offset edges, induced by pi-1 pi and pi pi+1, are connected by a circular arc of radius r, whose supporting circle is centered at pi. The angle that defines such a circular arc equals 180° - (pi-1, pi, pi+1); see Figure 22.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 P is not convex, its offset can be obtained by decomposing it to convex sub-polygons P1, Pm such that i=1mPi = P, computing the offset of each sub-polygon and finally calculating the union of these sub-offsets (see Figure 22.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 Br,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 pi is defined by an angle 540° - (pi-1, pi, pi+1); see Figure 22.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.

22.3.1   Approximating the Offset with a Guaranteed Error-Bound

Let us assume that we are given a counterclockwise-oriented polygon P = ( p0, , pn-1 ), whose vertices all have rational coordinates (i.e., for each vertex pi = (xi, yi) we have xi, yi , and we wish to compute its Minkowski sum with a disc of radius r, where r is also a rational number. The boundary of this sum is comprised of line segments and circular arcs, where:

Approximating an offset edge
Figure 22.5:  Approximating the offset edge o1 o2 induced by the polygon edge p1 p2 by two line segments o'1 q' and q' o'2.

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:

  1. Consider the example depicted in Figure 22.5, where the exact offset edge o1 o2 is obtained by shifting the polygon edge p1 p2 by a vector whose length is r that form an angle Φ with the x-axis. We select two points o'1 and o'2 with rational coordinates on the two circles centered at p1 and p2, respectively. These points are selected such that if we denote the angle that the vector pj oj forms with the x-axis by Φ'j (for j = 1, 2), we have Φ'1 < Φ< Φ'2.
  2. We construct two tangents to the two circles at o'1 and o'2, respectively. The tangent lines have rational coefficients.
  3. We compute the intersection point of the two tangents, denoted q'. The two line segments o'1 q' and q' o'2 approximate the original offset edge o1 o2.

The function approximated_offset_2 (P, r, eps) accepts a polygon P, an offset radius r and ε> 0. It constructs an approximation for the offset of P by the radius r using the procedure explained above. Furthermore, it is guaranteed that the approximation error, namely the distance of the point q' from o1 o2 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 P is not convex, its offset may not be simple and may contain holes, whose boundary is also comprised of line segments and circular arcs).

Offsetting a polygon
Figure 22.6:  The offset computation performed by the example programs Minkowski_sum_2/approx_offset.cpp and Minkowski_sum_2/exact_offset.cpp. The input polygon is shaded and the boundary of its offset is drawn in a thick black line.

The following example demonstrates the construction of an approximated offset of a non-convex polygon, as depicted in Figure 22.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;

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);
}

22.3.2   Computing the Exact Offset

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 (P, r, traits) computes the offset of the polygon P by a rational radius r 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 22.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;

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

Both functions approximated_offset_2() and offset_polygon_2() also have overloaded versions that accept a decomposition strategy and use the polygon-decomposition approach to compute (or approximate) the offset. These functions are less efficient than their counterparts that employ the convolution approach, and are only included in the package for the sake of completeness.

22.3.3   Computing Inner Offsets

An operation closely related to the offset computation, is obtaining the inner offset of a polygon, or insetting it by a given radius. The inset of a polygon P with radius r is the set of points iside P whose distance from the polygon boundary, denoted ∂P, is at least r - namely: { p P | dist(p, ∂P) r }. Note that the resulting set may not be connected in case P is a non-convex polygon that has some narrow components, and thus is characterized by a (possibly empty) set of polygons whose edges are line segments and circular arcs of radius r.

The offset can be computed using the convolution method if we traverse the polygon in a clockwise orientation (and not in counterclockwise orientation, as done in case of ofsetting a polygon). As in case of the offset functions, the Minkowski-sum package contains two functions for insetting a simple polygon:

The function approximated_inset_2 (P, r, eps, oi) accepts a polygon P, an inset radius r and ε> 0. It constructs an approximation for the inset of P by the radius r, with the approximation error bounded by ε. The function returns its output via the output iterator oi, whose value-type must be Gps_circle_segment_traits_2<Kernel>::Polygon_2 representing the polygons that approximates the inset polygon.

File: examples/Minkowski_sum_2/approx_inset.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;

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 ("tight.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 = 1;
  const double                           err_bound = 0.00001;
  std::list<Offset_polygon_2>            inset_polygons;
  std::list<Offset_polygon_2>::iterator  iit;
  CGAL::Timer                            timer;

  timer.start();
  approximated_inset_2 (P, radius, err_bound,
                        std::back_inserter (inset_polygons));
  timer.stop();

  std::cout << "The inset comprises " 
            << inset_polygons.size() << " polygon(s)." << std::endl;
  for (iit = inset_polygons.begin(); iit != inset_polygons.end(); ++iit)
  {
      std::cout << "    Polygon with "
                << iit->size() << " vertices." << std::endl;
  }
  std::cout << "Inset computation took "
            << timer.time() << " seconds." << std::endl;
  return (0);
}

Similarly, the function inset_polygon_2 (P, r, traits, oi) computes the exact inset of P with radius r, and returns its output via the given output iterator oi. The traits parameter should be an instance of an arrangement-traits class that is capable of handling conic arcs in an exact manner, whereas oi's value-type must be Gps_traits_2<ArrConicTraits>::Polygons_2.

File: examples/Minkowski_sum_2/exact_inset.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;

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;

int main ()
{
  // Open the input file.
  std::ifstream    in_file ("tight.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 inner offset of the polygon.
  Conic_traits_2                         traits;
  const Rational                         radius = 1;
  std::list<Offset_polygon_2>            inset_polygons;
  std::list<Offset_polygon_2>::iterator  iit;
  CGAL::Timer                            timer;

  timer.start();
  inset_polygon_2 (P, radius, traits,
                   std::back_inserter (inset_polygons));
  timer.stop();

  std::cout << "The inset comprises " 
            << inset_polygons.size() << " polygon(s)." << std::endl;
  for (iit = inset_polygons.begin(); iit != inset_polygons.end(); ++iit)
  {
      std::cout << "    Polygon with "
                << iit->size() << " vertices." << std::endl;
  }
  std::cout << "Inset computation took "
            << timer.time() << " seconds." << std::endl;
  return (0);
}

#endif

Unlike the offset functions, there are no overloaded versions of the inset functions that use convex polygon decomposition to compute insets, as this method cannot be easily generalized for inset computations.

In this context let us mention that there exist overloaded versions of the functions approximated_offset_2 (P, r, eps) and offset_polygon_2 (P, r, traits) that accept a polygon with holes P and computed its offset. This ofset is obtain by taking the outer offset of P's outer boundary, and computing the inner offsets of P's holes. The former polygon defines the output boundary of P Br, and the latter define the holes within the result.

Acknowledgements

Andreas Fabri and Laurent Rineau helped tracing and solving several bugs in the approximated offset function. They have also suggested a few algorithmic improvements that made their way into version 3.4, yielding a faster approximation scheme.


Footnotes

 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 v lies between two vectors u and w if we reach v strictly before reaching w if we move all three vectors to the origin and rotate u counterclockwise. Note that this also covers the case where u has the same direction as v.
 3  Informally speaking, the winding number of a point p 2 with respect to some planar curve γ is an integer number counting how many times does γ wind in a counterclockwise direction around p.
 4  As the disc is convex, it is guaranteed that the convolution curve is comprised of a single cycle.