CGAL 6.0.1 - 2D Straight Skeleton and Polygon Offsetting
Loading...
Searching...
No Matches
User Manual

Authors
Fernando Cacciola, Sébastien Loriot, and Mael Rouxel-Labbé

Introduction

Straight Skeletons

This package implements weighted straight skeletons for two-dimensional polygons with holes. An intuitive way to think of the construction of straight skeletons is to imagine that wavefronts (or grassfires) are spawned at each edge of the polygon, and are moving inward. As the fronts progress, they either contract or expand depending on the angles formed between polygon edges, and sometimes disappear. Under this transformation, polygon vertices move along the angular bisector of the lines subtending the edges, tracing a tree-like structure, the straight skeleton.

Figure 21.2 Construction of a straight skeleton: the wavefront interfaces define the straight skeleton.


Straight Skeletons, Medial Axis and Voronoi Diagrams

The straight skeleton of a polygon is similar to the medial axis and the Voronoi diagram of a polygon in the way it partitions it; however, unlike the medial axis and the Voronoi diagram, the bisectors are not equidistant to its defining edges but to the supporting lines of such edges. As a result, the bisectors of a straight skeleton might not be located in the center of the polygon and thus cannot be regarded as a proper medial axis in its geometrical meaning.

On the other hand, only reflex vertices (i.e., vertices whose internal angle \( > \pi\)) are responsible for deviations of the bisectors from its center location. Therefore, for convex polygons, the straight skeleton, the medial axis and the Voronoi diagram are exactly equivalent. Furthermore, if a non-convex polygon contains only vertices of low reflexivity, the straight skeleton bisectors will be placed nearly equidistant to their defining edges, producing a straight skeleton much alike a "proper" medial axis.

Uses of Straight Skeletons

This package implements construction of straight skeletons as well as two typical use cases of the straight skeleton: polygon offsetting and straight skeleton extrusion.

Definitions

We now define in detail what is the straight skeleton, starting from the input polygon.

2D Contours

A 2D contour is a closed sequence (a cycle) of three or more connected 2D oriented straight line segments called contour edges. The endpoints of contour edges are called vertices, and are shared by exactly two incident edges.

A contour partitions the plane in open regions that are bounded and unbounded. If the contour forms the boundary of a single open region that is a simply connected set, then this contour is said to be weakly simple. If the edges of the contour intersect only at their common vertices, the contour is said to be simple.

The orientation of a contour is given by the order of the vertices around the region they bound. It can be clockwise (CW) or counterclockwise (CCW). The bounded side of a contour edge is the side facing the bounded region of the contour. If the contour is oriented CCW, the bounded side of an edge is its left side.

2D Polygon with Holes

A 2D polygon is a contour.

A 2D polygon with holes is a contour having zero or more contours in its bounded regions. The former contour is called outer contour or outer boundary, and the latter contours are called inner contours, or holes. In this chapter, we require that there cannot be any intersection among any two contours. Furthermore, we require that the outer contour is CCW oriented, that holes are CW oriented, and that a hole cannot be in the bounded region of another hole.

The intersection of the bounded regions of the outer contour and the unbounded regions of each inner contour is the interior of the polygon with holes. A polygon with holes is said to be weakly simple if its interior is a simply connected set.

Throughout the rest of this chapter the term polygon will be used as a shortcut for polygon with holes.

Figure 21.3 Examples of weakly simple polygons: one with no holes and two edges coincident (left) and one with 2 holes (right).


Figure 21.4 Examples of non-simple polygons: one folding into itself, that is, non-planar (left), one with a vertex touching an edge (middle), and one with a hole crossing into the outside (right)


Straight Skeleton of a 2D Weakly Simple Polygon

Given a contour edge called the source edge, its offset edge at time t is an edge parallel with the same orientation such that the Euclidean distance between the lines supporting both edges is exactly t. An offset edge is always located to the bounded side of its source edge (which is an oriented straight line segment).

The straight skeleton is spawned by applying a continuous, inward offsetting to the contour edges. Under this transformation, polygon vertices move along the angular bisector of the lines subtending the edges, at a speed that depends on the angle between the two incident contour edges. So-called events happen when two moving vertices impact each other, or when a moving vertex collides with an offset edge. The different event types are detailed in the documentation of the class CGAL::Straight_skeleton_builder_2.

The straight skeleton is the set of segments traced out by the moving vertices during this process. It forms a unique partitioning of the polygon interior into straight skeleton regions corresponding to the monotone areas traced by the continuous inward offsetting of the contour edges. Each region corresponds to exactly 1 contour edge. These regions are bounded by the angular bisectors of the supporting lines of the contour edges, and each such region is itself a non-convex, weakly simple polygon.

Figure 21.5 Straight skeletons of various weakly simple polygons.


The main entry points for straight skeletons are the following functions:

Weighted Straight Skeletons

Weighted straight skeletons are a generalization of straight skeletons: contour edges are assigned a positive weight, which can be understood as assigning a speed to the wavefront spawned from the contour edge. Vertices no longer move along the angular bisector between two contour edges, but along a weighted bisector.

This CGAL package supports positive multiplicative weights: if the supporting line of a contour edge is described through the equation ax+by+c=0 then the supporting line of the offset edge at distance t is ax+by+c-t=0. With a multiplicative weight w, the equation becomes w(ax+by+c)-t=0. Therefore, a larger weight implies a faster moving front.

Figure 21.6 An unweighted straight skeleton (leftmost) and three randomly weighted straight skeletons.


A significant change from unweighted straight skeleton of polygon with holes is that the faces of a straight skeleton of a polygon with holes are no longer necessarily weakly simple polygons: a face can for example completely encompass a set of faces incident to a hole, giving it the topology of a ring. This property is essential, both to have a valid halfedge data structure (e.g., one cannot have a face with multiple boundaries), and to correctly extract the offset contours of the straight skeleton. Thus, a post processing step is applied to ensure the weakly simple property in all faces. This post processing adds so-called artificial bisector/vertices where required to cut the face into a weakly simple polygon. This is achieved by shooting a ray from the farthest (w.r.t. time) vertex for each boundary that does not contain the defining contour edge of its incident face.

Figure 21.7 A polygon with four holes (black) and its straight skeleton (red). Three artificial bisectors (green) are added in a post-processing step to recover the simply-connected property for the straight skeleton face (gray) of the rightmost vertical contour edge.


The main entry points for weighted straight skeletons are the following functions:

Straight Skeleton of a General Polygon

A straight skeleton can also be defined for a general multiply-connected, planar, directed straight-line graph by considering all the edges as embedded in an unbounded region [1]. The only difference is that in this case some faces will be only partially bounded.

The current version of this CGAL package can only construct the straight skeleton in the interior of a simple polygon with holes, that is it does not handle general polygonal figures in the plane.

Uses of Straight Skeletons

Polygon Offsetting

The straight skeleton is defined by the trace of the moving vertices, but one can also consider the state of the polygon at a fixed time t made up of the translated vertices. This is the offset polygon.

An offset polygon can have fewer, equal, or more sides as its source polygon. It can even be composed of multiple polygons. If the source polygon has no holes, then no offset polygon has any holes. If the source polygon has holes, any of the offset polygons can have holes itself, but it might as well have no holes at all (if the distance is sufficiently large). Each offset polygon has the same orientation as the source polygon.

The main entry points for polygon offsetting are the following functions:

Exterior Straight Skeletons and Exterior Offset Contours

This CGAL package can only construct the straight skeleton and offset contours in the interior of a polygon with holes. However, constructing exterior skeletons and exterior offsets is possible.

Say you have some polygon with one hole, and you wish to obtain some exterior offset contours. The interior region of a polygon with holes is connected while the exterior region is not: there is an unbounded region outside the outer contour, and one bounded region inside each hole. To construct an offset contour you need to construct a straight skeleton. Thus, to construct exterior offset contours for a polygon with holes, you need to construct, separately, the exterior skeleton of the outer contour and the interior skeleton of each hole.

Constructing the interior skeleton of a hole is directly supported by this CGAL package; you simply need to input the hole's vertices in reversed order as if it were an outer contour as to treat them as simple polygons without holes. The offset contour(s) should also be reversed.

Constructing the exterior skeleton of the outer contour is achieved by means of the following trick: consider the outer contour as a hole of a larger rectangle (call it frame). If the frame is sufficiently separated from the contour, the resulting skeleton will be equivalent to a real exterior skeleton once the offset contour corresponding to the outer frame is removed, which is performed automatically by some of the functions of this package.

Figure 21.8 Exterior skeleton obtained using a frame (left) and 2 sample exterior offset contours (right)


Advanced

It is necessary to place the frame sufficiently far away from the contour. If it is not, it could occur that the outward offset contour collides and merges with the inward offset frame, resulting in 1 instead of 2 offset contours. However, the proper separation between the contour and the frame is not directly given by the offset distance at which you want the offset contour. That distance must be at least the desired offset plus the largest euclidean distance between an offset vertex and its original. This CGAL packages provides a helper function to compute the required separation: compute_outer_frame_margin().

For convenience, the following functions are provided:

These also exists in weighted versions.

Straight Skeleton Extrusion

Perhaps the first (historically) use-case of straight skeletons: given a polygonal roof, the straight skeleton directly gives the layout of each tent. If each skeleton edge is lifted from the plane a height equal to its offset distance, the resulting roof is "correct" in that water will always fall down to the contour edges (the roof's border), regardless of where it falls on the roof. Laycock and Day [4] give an algorithm for roof design based on the straight skeleton.

This CGAL package implements skeleton extrusion for polygons with holes, with support for positive multiplicative weights. The output is a closed, combinatorially 2-manifold surface triangle mesh.

Figure 21.9 Input polygon (left), weighted skeleton with colored faces (middle), and extruded skeleton (right).


The main entry point for straight skeleton extrusion is the function CGAL::extrude_skeleton(), whose API also enables passing angles for each tent, which are then converted internally to edge weights.

The extrusion can be performed to a maximum height, as shown in the figure below.

Figure 21.10 Input polygon with three holes and its weighted straight skeleton (left), and two extrusions of the skeleton with different maximum heights (middle and right).


Other Uses of Straight Skeletons

Just like medial axes, straight skeletons can also be used for 2D shape description and matching. Essentially, all the applications of image-based skeletonization (for which there is a vast literature) are also direct applications of the straight skeleton, especially since skeleton edges are simply straight line segments.

Consider the subgraph formed only by inner bisectors (that is, only the skeleton edges which are not incident upon a contour vertex). Call this subgraph a skeleton axis. Each node in the skeleton axis whose degree is \( >=3\) roots more than one skeleton tree. Each skeleton tree roughly corresponds to a region in the input topologically equivalent to a rectangle; that is, without branches. For example, a simple character "H" would contain 2 higher degree nodes separating the skeleton axis in 5 trees; while the character "@" would contain just 1 higher degree node separating the skeleton axis in 2 curly trees.

Since a skeleton edge is a straight line, each branch in a skeleton tree is a polyline. Thus, the path-length of the tree can be directly computed. Furthermore, the polyline for a particular tree can be interpolated to obtain curve-related information.

Pruning each skeleton tree cutting off branches whose length is below some threshold; or smoothing a given branch, can be used to reconstruct the polygon without undesired details, or fit into a particular canonical shape.

Each skeleton edge in a skeleton branch is associated with 2 contour edges which are facing each other. If the polygon has a bottleneck (it almost touches itself), a search in the skeleton graph measuring the distance between each pair of contour edges will reveal the location of the bottleneck, allowing you to cut the shape in two. Likewise, if two shapes are too close to each other along some part of their boundaries (a near contact zone), a similar search in an exterior skeleton of the two shapes at once would reveal the parts of near contact, allowing you to stitch the shapes. These cut and stitch operations can be directly executed in the straight skeleton itself instead of the input polygon (because the straight skeleton contains a graph of the connected contour edges).

Implementation Details

The reference manual of this CGAL package contains many details about the straight skeleton construction algorithm that has been implemented. In particular, the classes and CGAL::Straight_skeleton_2 and CGAL::Straight_skeleton_builder_2 describe the algorithm in depth.

Examples

Create a Straight Skeleton

The straight skeleton data structure is implemented in the class CGAL::Straight_skeleton_2.

The simplest way to construct a straight skeleton is via the free functions CGAL::create_interior_straight_skeleton_2() and CGAL::create_exterior_straight_skeleton_2(), as shown in the following example:


File Straight_skeleton_2/Create_straight_skeleton_2.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_2.h>
#include <CGAL/create_straight_skeleton_2.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <CGAL/draw_straight_skeleton_2.h>
#include <cassert>
typedef K::Point_2 Point ;
typedef CGAL::Polygon_2<K> Polygon_2 ;
typedef std::shared_ptr<Ss> SsPtr ;
int main()
{
Polygon_2 poly ;
poly.push_back( Point(-1,-1) ) ;
poly.push_back( Point(0,-12) ) ;
poly.push_back( Point(1,-1) ) ;
poly.push_back( Point(12,0) ) ;
poly.push_back( Point(1,1) ) ;
poly.push_back( Point(0,12) ) ;
poly.push_back( Point(-1,1) ) ;
poly.push_back( Point(-12,0) ) ;
assert(poly.is_counterclockwise_oriented());
// You can pass the polygon via an iterator pair
SsPtr iss = CGAL::create_interior_straight_skeleton_2(poly.vertices_begin(), poly.vertices_end());
CGAL::Straight_skeletons_2::IO::print_straight_skeleton(*iss);
CGAL::draw(*iss);
// Or you can pass the polygon directly, as below.
// To create an exterior straight skeleton you need to specify a maximum offset.
double lMaxOffset = 5 ;
SsPtr oss = CGAL::create_exterior_straight_skeleton_2(lMaxOffset, poly);
CGAL::Straight_skeletons_2::IO::print_straight_skeleton(*oss);
CGAL::draw(*oss);
return EXIT_SUCCESS;
}
The class Straight_skeleton_2 provides a model for the StraightSkeleton_2 concept which is the class ...
Definition: Straight_skeleton_2.h:178
void draw(const PS2 &ps2, const GSOptions &gso)
std::shared_ptr< Straight_skeleton_2< SsK > > create_exterior_straight_skeleton_2(FT max_offset, PointIterator vertices_begin, PointIterator vertices_end, SsK k=CGAL::Exact_predicates_inexact_constructions_kernel())
creates a straight skeleton in the exterior of a 2D polygon with holes.
std::shared_ptr< Straight_skeleton_2< SsK > > create_interior_straight_skeleton_2(PointIterator outer_contour_vertices_begin, PointIterator outer_contour_vertices_end, HoleIterator holes_begin, HoleIterator holes_end, SsK k=CGAL::Exact_predicates_inexact_constructions_kernel())
creates a straight skeleton in the interior of the 2D polygon with holes.

The input to these functions is the polygon, which can be given as an iterator pair or directly as a CGAL::Polygon_2 object. In the case of the exterior skeleton, a maximum offset must be specified as well (see Section Straight_skeleton_2ExteriorSkeletonsandExterior for details on this max offset parameter).

Create a Straight Skeleton from a Polygon With Holes

If CGAL::Polygon_with_holes_2 is used, you can pass an instance of it directly to the function creating the interior skeleton, as shown below. Notice that a different header must be included in this case.


File Straight_skeleton_2/Create_straight_skeleton_from_polygon_with_holes_2.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_with_holes_2.h>
#include <CGAL/create_straight_skeleton_from_polygon_with_holes_2.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <memory>
#include <cassert>
typedef K::Point_2 Point ;
typedef CGAL::Polygon_2<K> Polygon_2 ;
typedef CGAL::Polygon_with_holes_2<K> Polygon_with_holes ;
typedef std::shared_ptr<Ss> SsPtr ;
int main()
{
Polygon_2 outer ;
outer.push_back( Point(-1,-1) ) ;
outer.push_back( Point(0,-12) ) ;
outer.push_back( Point(1,-1) ) ;
outer.push_back( Point(12,0) ) ;
outer.push_back( Point(1,1) ) ;
outer.push_back( Point(0,12) ) ;
outer.push_back( Point(-1,1) ) ;
outer.push_back( Point(-12,0) ) ;
Polygon_2 hole ;
hole.push_back( Point(-1,0) ) ;
hole.push_back( Point(0,1 ) ) ;
hole.push_back( Point(1,0 ) ) ;
hole.push_back( Point(0,-1) ) ;
assert(outer.is_counterclockwise_oriented());
assert(hole.is_clockwise_oriented());
Polygon_with_holes poly( outer ) ;
poly.add_hole( hole ) ;
CGAL::Straight_skeletons_2::IO::print_straight_skeleton(*iss);
return EXIT_SUCCESS;
}

Create Offset Polygons from a Straight Skeleton

If you already have a straight skeleton object, the simpler way to generate offset polygons is to call CGAL::create_offset_polygons_2() as shown in the next example, passing the desired offset and the straight skeleton. You can reuse the same skeleton to generate offsets at a different distance, which is recommended because producing the straight skeleton is much slower than generating offset polygons.


File Straight_skeleton_2/Create_offset_polygons_2.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_2.h>
#include <CGAL/create_offset_polygons_2.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <memory>
#include <vector>
#include <cassert>
typedef K::Point_2 Point ;
typedef CGAL::Polygon_2<K> Polygon_2 ;
typedef std::shared_ptr<Polygon_2> PolygonPtr ;
typedef std::shared_ptr<Ss> SsPtr ;
typedef std::vector<PolygonPtr> PolygonPtrVector ;
int main()
{
Polygon_2 poly ;
poly.push_back( Point(-1,-1) ) ;
poly.push_back( Point(0,-12) ) ;
poly.push_back( Point(1,-1) ) ;
poly.push_back( Point(12,0) ) ;
poly.push_back( Point(1,1) ) ;
poly.push_back( Point(0,12) ) ;
poly.push_back( Point(-1,1) ) ;
poly.push_back( Point(-12,0) ) ;
assert(poly.is_counterclockwise_oriented());
double lOffset = 1 ;
PolygonPtrVector offset_polygons = CGAL::create_offset_polygons_2<Polygon_2>(lOffset,*ss);
CGAL::Straight_skeletons_2::IO::print_polygons(offset_polygons);
return EXIT_SUCCESS;
}

Create Offset Polygons from a Polygon (With or Without Holes)

If you need offset polygons at a single distance, you can hide away the construction of the straight skeleton by calling directly the functions CGAL::create_interior_skeleton_and_offset_polygons_2() and CGAL::create_exterior_skeleton_and_offset_polygons_2() as shown in the following examples:


File Straight_skeleton_2/Create_skeleton_and_offset_polygons_2.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_2.h>
#include <CGAL/create_offset_polygons_2.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <memory>
#include <cassert>
#include <vector>
typedef K::FT FT ;
typedef K::Point_2 Point ;
typedef CGAL::Polygon_2<K> Polygon_2 ;
typedef std::shared_ptr<Polygon_2> PolygonPtr ;
typedef std::shared_ptr<Ss> SsPtr ;
typedef std::vector<PolygonPtr> PolygonPtrVector ;
int main()
{
Polygon_2 poly ;
poly.push_back( Point(-1,-1) ) ;
poly.push_back( Point(0,-12) ) ;
poly.push_back( Point(1,-1) ) ;
poly.push_back( Point(12,0) ) ;
poly.push_back( Point(1,1) ) ;
poly.push_back( Point(0,12) ) ;
poly.push_back( Point(-1,1) ) ;
poly.push_back( Point(-12,0) ) ;
assert(poly.is_counterclockwise_oriented());
FT lOffset = 1 ;
PolygonPtrVector inner_offset_polygons = CGAL::create_interior_skeleton_and_offset_polygons_2(lOffset,poly);
PolygonPtrVector outer_offset_polygons = CGAL::create_exterior_skeleton_and_offset_polygons_2(lOffset,poly);
CGAL::Straight_skeletons_2::IO::print_polygons(inner_offset_polygons);
CGAL::Straight_skeletons_2::IO::print_polygons(outer_offset_polygons);
return EXIT_SUCCESS;
}
std::vector< std::shared_ptr< OfKPolygon > > create_interior_skeleton_and_offset_polygons_2(FT offset, const InKPolygon &outer_boundary, HoleIterator holes_begin, HoleIterator holes_end, OfK ofk=CGAL::Exact_predicates_inexact_constructions_kernel, SsK ssk=CGAL::Exact_predicates_inexact_constructions_kernel())
returns a container with all the inner offset polygons at distance offset of a 2D polygon with holes.
std::vector< std::shared_ptr< OfKPolygon > > create_exterior_skeleton_and_offset_polygons_2(FT offset, const InKPolygon &poly, OfK ofk=Exact_predicates_inexact_constructions_kernel(), SsK ssk=Exact_predicates_inexact_constructions_kernel())
returns a container with all the outer offset polygons at distance offset of the 2D polygon poly.

... and using a Polygon_with_holes_2 directly when available:


File Straight_skeleton_2/Create_saop_from_polygon_with_holes_2.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_with_holes_2.h>
#include <CGAL/create_offset_polygons_from_polygon_with_holes_2.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <memory>
#include <cassert>
#include <vector>
typedef K::Point_2 Point ;
typedef CGAL::Polygon_2<K> Polygon_2 ;
typedef CGAL::Polygon_with_holes_2<K> Polygon_with_holes ;
typedef std::shared_ptr<Polygon_2> PolygonPtr ;
typedef std::shared_ptr<Ss> SsPtr ;
typedef std::vector<PolygonPtr> PolygonPtrVector ;
int main()
{
Polygon_2 outer ;
outer.push_back( Point(-1,-1) ) ;
outer.push_back( Point(0,-12) ) ;
outer.push_back( Point(1,-1) ) ;
outer.push_back( Point(12,0) ) ;
outer.push_back( Point(1,1) ) ;
outer.push_back( Point(0,12) ) ;
outer.push_back( Point(-1,1) ) ;
outer.push_back( Point(-12,0) ) ;
Polygon_2 hole ;
hole.push_back( Point(-1,0) ) ;
hole.push_back( Point(0,1 ) ) ;
hole.push_back( Point(1,0 ) ) ;
hole.push_back( Point(0,-1) ) ;
assert(outer.is_counterclockwise_oriented());
assert(hole.is_clockwise_oriented());
Polygon_with_holes poly( outer ) ;
poly.add_hole( hole ) ;
double lOffset = 0.2 ;
PolygonPtrVector offset_polygons = CGAL::create_interior_skeleton_and_offset_polygons_2(lOffset,poly);
CGAL::Straight_skeletons_2::IO::print_polygons(offset_polygons);
return EXIT_SUCCESS;
}

If the input polygon has holes, there can be holes in the offset polygons. However, the polygons generated by all the offsetting functions shown before do not have any parent-hole relationship computed; that is, they just instances of Polygon_2 instead of Polygon_with_holes_2. If Polygon_with_holes_2 are available and you need the offsetting to produce them, you can call the function CGAL::arrange_offset_polygons_2() passing the result of any of the offsetting functions described so far. That function arranges the offset polygons detecting and distributing holes within parents. As a shortcut, you can use the function CGAL::create_interior_skeleton_and_offset_polygons_with_holes_2() as shown below:


File Straight_skeleton_2/Create_skeleton_and_offset_polygons_with_holes_2.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_with_holes_2.h>
#include <CGAL/create_offset_polygons_from_polygon_with_holes_2.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <memory>
#include <vector>
#include <cassert>
typedef K::Point_2 Point ;
typedef CGAL::Polygon_2<K> Polygon_2 ;
typedef CGAL::Polygon_with_holes_2<K> PolygonWithHoles ;
typedef std::shared_ptr<PolygonWithHoles> PolygonWithHolesPtr ;
typedef std::vector<PolygonWithHolesPtr> PolygonWithHolesPtrVector;
int main()
{
Polygon_2 outer ;
outer.push_back( Point( 0.0, 0.0) ) ;
outer.push_back( Point(10.0, 0.0) ) ;
outer.push_back( Point(10.0, 4.5) ) ;
outer.push_back( Point(12.0, 4.5) ) ;
outer.push_back( Point(12.0, 2.0) ) ;
outer.push_back( Point(16.0, 2.0) ) ;
outer.push_back( Point(16.0, 8.0) ) ;
outer.push_back( Point(12.0, 8.0) ) ;
outer.push_back( Point(12.0, 5.5) ) ;
outer.push_back( Point(10.0, 5.5) ) ;
outer.push_back( Point(10.0,10.0) ) ;
outer.push_back( Point( 0.0,10.0) ) ;
Polygon_2 hole ;
hole.push_back( Point(3.0,3.0) ) ;
hole.push_back( Point(3.0,7.0) ) ;
hole.push_back( Point(7.0,7.0) ) ;
hole.push_back( Point(7.0,3.0) ) ;
assert(outer.is_counterclockwise_oriented());
assert(hole.is_clockwise_oriented());
PolygonWithHoles poly( outer ) ;
poly.add_hole( hole ) ;
double lOffset = 1 ;
PolygonWithHolesPtrVector offset_poly_with_holes = CGAL::create_interior_skeleton_and_offset_polygons_with_holes_2(lOffset,poly);
CGAL::Straight_skeletons_2::IO::print_polygons_with_holes(offset_poly_with_holes);
return EXIT_SUCCESS;
}
std::vector< std::shared_ptr< OfKPolygon > > create_interior_skeleton_and_offset_polygons_with_holes_2(FT offset, const InKPolygon &poly_with_holes, OfK ofk=CGAL::Exact_predicates_inexact_constructions_kernel, SsK ssk=CGAL::Exact_predicates_inexact_constructions_kernel)
returns a container with all the inner offset polygons with holes at distance offset of the 2D polygo...
Advanced

Consider an input polygon with parallel edges separated a distance \( 2*t\). If you produce an offset polygon at distance \( t\), these parallel edges will just collapse each other and vanish from the result, keeping the output as a simple polygon, just like the input. However, if you request an offset polygon at a distance \( t-epsilon\), the result will still be a simple polygon but with edges that are so close to each other that will almost intersect. If a kernel with exact constructions is used, the offsetting algorithm can guarantee that the output contains only simple polygons. However, if inexact constructions are used the roundoff in the coordinates of the output points will cause parallel edges that almost collapse-but not so-to become really collinear or even cross each other.

Thus, it is necessary to use a kernel with exact constructions if offset polygons must be simple, yet computing a straight skeleton using that kernel is very slow, much more than computing the offset polygons. To help with this, it is possible to construct the straight skeleton using the recommended kernel Exact_predicates_inexact_constructions_kernel, then convert the skeleton to a different kernel via the function CGAL::convert_straight_skeleton_2() and input the converted skeleton to the offsetting functions.

All the offsetting functions that take polygons as input (and create the straight skeleton under the hood) apply that optimization automatically: that is, the output polygons are defined over the same kernel of the input polygons, whatever that is, yet the straight skeleton is constructed with the faster recommended kernel and converted if necessary.

Notice how some of the examples above use Exact_predicates_exact_constructions_kernel. In all cases, the straight skeleton is constructed using Exact_predicates_inexact_constructions_kernel.

Extrude the Skeleton of a Polygon with Holes

The following example shows how to extrude the weighted skeleton of a polygon with holes:


File Straight_skeleton_2/extrude_skeleton.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
//#include <CGAL/Exact_predicates_exact_constructions_kernel_with_sqrt.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/draw_straight_skeleton_2.h>
#include <CGAL/draw_surface_mesh.h>
#include "CGAL/input_helpers.h" // polygon reading, random polygon with weights generation
#include <CGAL/extrude_skeleton.h>
#include <CGAL/IO/polygon_mesh_io.h>
#include <CGAL/Polygon_2.h>
#include <CGAL/Real_timer.h>
#include <CGAL/Random.h>
#include <iostream>
#include <unordered_map>
#include <vector>
#include <memory>
namespace SS = CGAL::CGAL_SS_i;
namespace PMP = CGAL::Polygon_mesh_processing;
// Kernel choice:
// EPICK: Robust and fast
// EPECK_with_sqrt: Exact and slow
// EPECK: More robust, and less slow than EPECK_with_sqrt
// using K = CGAL::Exact_predicates_exact_constructions_kernel;
// using K = CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt;
using FT = K::FT;
using Point_2 = K::Point_2;
using Segment_2 = K::Segment_2;
using Line_2 = K::Line_2;
using Point_3 = K::Point_3;
using Vector_3 = K::Vector_3;
using Polygon_2 = CGAL::Polygon_2<K>;
using Polygon_with_holes_2 = CGAL::Polygon_with_holes_2<K>;
using Straight_skeleton_2 = CGAL::Straight_skeleton_2<K>;
using Straight_skeleton_2_ptr = std::shared_ptr<Straight_skeleton_2>;
using Mesh = CGAL::Surface_mesh<Point_3>;
int main(int argc, char** argv)
{
std::cout.precision(17);
std::cerr.precision(17);
const int argc_check = argc - 1;
char* poly_filename = nullptr;
char* speeds_filename = nullptr;
FT height = FT{(std::numeric_limits<double>::max)()};
bool use_angles = false; // whether the input is SLS edge weights, or taper angles
bool flip_weights = false; // takes the opposite for weights, and the complement for angles
// below is only used for random weight generation
double min_weight = 1., max_weight = 10.;
int seed = CGAL::get_default_random().get_seed();
for(int i = 1; i < argc; ++i)
{
if(!strcmp("-h", argv[i]) || !strcmp("--help", argv[i]) || !strcmp("-?", argv[i]))
{
std::cout << "Usage: " << argv[0] << "[options].\n"
"Options:\n"
" -i <input_filename>: input polygon filename.\n"
" -t <value>: height. Must be non-zero.\n"
" -a <angles_filename>: angles. Format: one angle per line, a space to separate borders.\n"
" -w <weights_filename>: weights. Format: one weight per line, a space to separate borders.\n"
" Note: -w and -a are exclusive.\n"
<< std::endl;
return EXIT_FAILURE;
} else if(!strcmp("-i", argv[i]) && i < argc_check) {
poly_filename = argv[++i];
} else if(!strcmp("-w", argv[i])) {
if(speeds_filename != nullptr)
{
std::cerr << "Error: -w and -a are exclusive." << std::endl;
return EXIT_FAILURE;
}
speeds_filename = argv[++i];
use_angles = false;
} else if(!strcmp("-a", argv[i])) {
if(speeds_filename != nullptr)
{
std::cerr << "Error: -w and -a are exclusive." << std::endl;
return EXIT_FAILURE;
}
speeds_filename = argv[++i];
use_angles = true;
} else if(!strcmp("-t", argv[i]) && i < argc_check) {
height = std::stod(argv[++i]);
} else if(!strcmp("-mw", argv[i]) && i < argc_check) {
min_weight = std::stod(argv[++i]);
} else if(!strcmp("-Mw", argv[i]) && i < argc_check) {
max_weight = std::stod(argv[++i]);
} else if(!strcmp("-f", argv[i]) && i < argc) {
flip_weights = true;
} else if(!strcmp("-s", argv[i]) && i < argc_check) {
seed = std::stoi(argv[++i]);
}
}
if(CGAL::is_zero(height))
{
std::cerr << "Error: height must be non-zero" << std::endl;
return EXIT_FAILURE;
}
CGAL::Real_timer timer;
timer.start();
Polygon_with_holes_2 pwh;
if(poly_filename == nullptr)
{
pwh = generate_random_polygon<Polygon_with_holes_2>(seed);
}
else if(!read_input_polygon(poly_filename, pwh) || pwh.outer_boundary().is_empty())
{
std::cerr << "Error: failure during polygon read" << std::endl;
return EXIT_FAILURE;
}
#ifdef CGAL_SLS_OUTPUT_FILES
std::ofstream out_poly("input.dat");
out_poly.precision(17);
out_poly << pwh;
out_poly.close();
#endif
// read segment speeds (angles or weights)
std::vector<std::vector<FT> > speeds;
if(speeds_filename == nullptr)
generate_random_weights(pwh, min_weight, max_weight, seed, speeds);
else
read_segment_speeds(speeds_filename, speeds);
if(flip_weights)
{
if(use_angles)
{
for(auto& contour_speeds : speeds)
for(FT& a : contour_speeds)
a = 180 - a;
}
else
{
for(auto& contour_speeds : speeds)
for(FT& w : contour_speeds)
w = -w;
}
}
timer.stop();
std::cout << "Reading input(s) took " << timer.time() << " s." << std::endl;
// End of I/O, do some slope preprocessing and check the validity of the input(s)
// -----------------------------------------------------------------------------------------------
timer.reset();
timer.start();
Mesh sm;
if(use_angles)
CGAL::extrude_skeleton(pwh, sm, CGAL::parameters::angles(speeds).maximum_height(height));
else
CGAL::extrude_skeleton(pwh, sm, CGAL::parameters::weights(speeds).maximum_height(height));
timer.stop();
std::cout << "Extrusion computation took " << timer.time() << " s." << std::endl;
CGAL::IO::write_polygon_mesh("extruded_skeleton.off", sm, CGAL::parameters::stream_precision(17));
return EXIT_SUCCESS;
}
result_type is_zero(const NT &x)
bool write_polygon_mesh(const std::string &fname, Graph &g, const NamedParameters &np=parameters::default_values())
bool extrude_skeleton(const Polygon &polygon, PolygonMesh &out, const NamedParameters &np=parameters::default_values())
constructs the straight skeleton-based extrusion of a polygon with holes.

Low Level API

All the high level functions described above are just wrappers around the low-level API described here. This low level API is richer and provides options and configurations not covered by any of those functions.

The straight skeleton construction algorithm is encapsulated in the class CGAL::Straight_skeleton_builder_2, which is parameterized by a geometric traits (CGAL::Straight_skeleton_builder_traits_2) and the Straight Skeleton class.

The offset contours construction algorithm is encapsulated in the class CGAL::Polygon_offset_builder_2, which is parameterized by the Straight Skeleton class, a geometric traits (Polygon_offset_builder_traits_2<Kernel>) and a container type where the resulting offset polygons are generated.

To construct the straight skeleton of a polygon the user must:

  1. Instantiate the straight skeleton builder.
  2. Enter one contour at a time, starting from the outer contour, via the method Straight_skeleton_builder_2::enter_contour(). The input polygon must be weakly simple and counter-clockwise oriented (see the definitions at the beginning of this chapter). Collinear edges are allowed. The insertion order of each hole is unimportant but the outer contour must be entered first.
  3. Call Straight_skeleton_builder_2::construct_skeleton() once all the contours have been entered. You cannot enter another contour once the skeleton has been constructed.

To construct a set of inward offset contours the user must:

  1. Construct the straight skeleton of the source polygon.
  2. Instantiate the polygon offset builder passing in the straight skeleton as a parameter.
  3. Call Polygon_offset_builder_2::construct_offset_contours() passing the desired offset distance and an output iterator that can store a std::shared_ptr of Container instances into a resulting sequence (typically, a back insertion iterator)

Each element in the resulting sequence is an offset contour, given by a std::shared_ptr holding a dynamically allocated instance of the Container type. Such a container can be any model of the SequenceContainer concept, for example, a CGAL::Polygon_2, or just a std::vector of 2D points.

The resulting sequence of offset contours can contain both outer and inner contours. Each offset hole (inner offset contour) would logically belong in the interior of some of the outer offset contours. However, this algorithm returns a sequence of contours in arbitrary order and there is no indication whatsoever of the parental relationship between inner and outer contours.

On the other hand, each outer contour is counter-clockwise oriented while each hole is clockwise-oriented. And since offset contours do form simple polygons, it is guaranteed that no hole will be inside another hole, no outer contour will be inside any other contour, and each hole will be inside exactly 1 outer contour.

Parental relationships are not automatically reconstructed by this algorithm because this relation is not directly given by the input polygon and must be done as a post-processing step. The function CGAL::arrange_offset_polygons_2() can be used to do that efficiently.

A user can reconstruct the parental relationships as a post processing operation by testing each inner contour (which is identified by being clockwise) against each outer contour (identified as being counter-clockwise) for insideness.

This algorithm requires exact predicates but not exact constructions. Therefore, the Exact_predicates_inexact_constructions_kernel should be used.


File Straight_skeleton_2/Low_level_API.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polygon_2.h>
#include <CGAL/Polygon_2_algorithms.h>
#include <CGAL/Straight_skeleton_builder_2.h>
#include <CGAL/Polygon_offset_builder_2.h>
#include <CGAL/compute_outer_frame_margin.h>
#include <CGAL/Straight_skeleton_2/IO/print.h>
#include <memory>
#include <vector>
#include <cassert>
//
// This example illustrates how to use the CGAL Straight Skeleton package
// to construct an offset contour on the outside of a polygon
//
// This is the recommended kernel
typedef Kernel::Point_2 Point_2;
typedef CGAL::Polygon_2<Kernel> Contour;
typedef std::shared_ptr<Contour> ContourPtr;
typedef std::vector<ContourPtr> ContourSequence ;
typedef Ss::Halfedge_iterator Halfedge_iterator;
typedef Ss::Halfedge_handle Halfedge_handle;
typedef Ss::Vertex_handle Vertex_handle;
int main()
{
// A start-shaped polygon, oriented counter-clockwise as required for outer contours.
Point_2 pts[] = { Point_2(-1,-1)
, Point_2(0,-12)
, Point_2(1,-1)
, Point_2(12,0)
, Point_2(1,1)
, Point_2(0,12)
, Point_2(-1,1)
, Point_2(-12,0)
} ;
std::vector<Point_2> star(pts,pts+8);
// We want an offset contour in the outside.
// Since the package doesn't support that operation directly, we use the following trick:
// (1) Place the polygon as a hole of a big outer frame.
// (2) Construct the skeleton on the interior of that frame (with the polygon as a hole)
// (3) Construct the offset contours
// (4) Identify the offset contour that corresponds to the frame and remove it from the result
double offset = 3 ; // The offset distance
// First we need to determine the proper separation between the polygon and the frame.
// We use this helper function provided in the package.
std::optional<double> margin = CGAL::compute_outer_frame_margin(star.begin(),star.end(),offset);
// Proceed only if the margin was computed (an extremely sharp corner might cause overflow)
if ( margin )
{
// Get the bbox of the polygon
CGAL::Bbox_2 bbox = CGAL::bbox_2(star.begin(),star.end());
// Compute the boundaries of the frame
double fxmin = bbox.xmin() - *margin ;
double fxmax = bbox.xmax() + *margin ;
double fymin = bbox.ymin() - *margin ;
double fymax = bbox.ymax() + *margin ;
// Create the rectangular frame
Point_2 frame[4]= { Point_2(fxmin,fymin)
, Point_2(fxmax,fymin)
, Point_2(fxmax,fymax)
, Point_2(fxmin,fymax)
} ;
// Instantiate the skeleton builder
SsBuilder ssb ;
// Enter the frame
ssb.enter_contour(frame,frame+4);
// Enter the polygon as a hole of the frame (NOTE: as it is a hole we insert it in the opposite orientation)
ssb.enter_contour(star.rbegin(),star.rend());
// Construct the skeleton
std::shared_ptr<Ss> ss = ssb.construct_skeleton();
// Proceed only if the skeleton was correctly constructed.
if ( ss )
{
CGAL::Straight_skeletons_2::IO::print_straight_skeleton(*ss);
// Instantiate the container of offset contours
ContourSequence offset_contours ;
// Instantiate the offset builder with the skeleton
OffsetBuilder ob(*ss);
// Obtain the offset contours
ob.construct_offset_contours(offset, std::back_inserter(offset_contours));
// Locate the offset contour that corresponds to the frame
// That must be the outmost offset contour, which in turn must be the one
// with the largetst unsigned area.
ContourSequence::iterator f = offset_contours.end();
double lLargestArea = 0.0 ;
for (ContourSequence::iterator i = offset_contours.begin(); i != offset_contours.end(); ++ i )
{
double lArea = CGAL_NTS abs( (*i)->area() ) ; //Take abs() as Polygon_2::area() is signed.
if ( lArea > lLargestArea )
{
f = i ;
lLargestArea = lArea ;
}
}
// Remove the offset contour that corresponds to the frame.
offset_contours.erase(f);
CGAL::Straight_skeletons_2::IO::print_polygons(offset_contours);
}
}
return 0;
}
double ymax() const
double xmin() const
double ymin() const
double xmax() const
The class Polygon_offset_builder_2 encapsulates the construction of inward offset contours of a 2D si...
Definition: Polygon_offset_builder_2.h:20
Definition: Polygon_offset_builder_traits_2.h:17
The class Straight_skeleton_builder_2 encapsulates the construction of the 2D straight skeleton in th...
Definition: Straight_skeleton_builder_2.h:118
Definition: Straight_skeleton_builder_traits_2.h:17
NT abs(const NT &x)
Orientation orientation_2(ForwardIterator first, ForwardIterator last, const Traits &traits)
std::optional< typename Traits::FT > compute_outer_frame_margin(InputIterator first, InputIterator beyond, typename Traits::FT offset, const Traits &traits=Default_traits)
computes the separation required between a polygon and the outer frame used to obtain an exterior ske...
const CGAL::Orientation COUNTERCLOCKWISE

Implementation History

Fernando Cacciola created the first version of this package, providing unweighted straight skeletons and polygon offsetting.

Sébastien Loriot and Mael Rouxel-Labbé implemented further robustification techniques, as well as weighted straight skeletons, and skeleton extrusion.