Chapter 31
Envelopes of Curves in 2D

Ron Wein

Table of Contents

31.1 Introduction
31.2 The Envelope Diagram
31.3 Examples
   31.3.1   Example for Envelope of Line Segments
   31.3.2   Example for Computing the Convex Hull with Envelopes
   31.3.3   Example for Envelope of Non-Linear Curves

Lower Envelope
Figure 31.1:  The lower envelope of a set of line segments and hyperbolic arc.

31.1   Introduction

A continuous curve C in 2 is called x-monotone, if every vertical line intersects it at a single point at most. For example, the circle x2 + y2 = 1 is not xy-monotone as the vertical line x = 0 intersects it at (0, -1) and at (0, 1); however, it is possible to split the circle into an upper part and a lower part, such that both of these parts are x-monotone. We consider vertical segments as weakly x-monotone, to properly handle inputs that contain such vertical curves.

An x-monotone curve can be represented as a univariate function y = C(x), defined over some continuous range RC . Given a set C = { C1, C2, , Cn } of x-monotone curves, their lower envelope is defined as the point-wise minimum of all curves. Namely, the lower envelope of the set C can be defined as the following function:

LC (x) = min 1 k nCk (x) ,
where we define Ck(x) = Ck(x) for x RCk, and Ck(x) = ∞ otherwise.

Similarly, the upper envelope of C is the point-wise maximum of the x-monotone curves in the set:

UC (x) = max 1 k nCk (x) ,
where in this case Ck(x) = -∞ for x RCk.

Given a set of x-monotone curves C, the minimization diagram of C is a subdivision of the x-axis into cells, such that the identity of the curves that induce the lower envelope over a specific cell of the subdivision (an edge or a vertex) is the same. In non-degenerate situations, an edge - which represents a continuous interval on the x-axis - is induced by a single curve (or by no curves at all, if there are no x-monotone curves defined over the interval), and a vertex is induced by a single curve and corresponds to one of its endpoints, or by two curves and corresponds to their intersection point. The maximization diagram is symmetrically defined for upper envelopes. In the rest of this chapter, we refer to both these diagrams as envelope diagrams.

Lower and upper envelopes can be efficiently computed using a divide-and-conquer approach. First, note that the envelope diagram for a single x-monotone curve Ck is trivial to compute: we project the boundary of its range of definition RCk onto the x-axis and label the features it induces accordingly. Given a set D of (non necessarily x-monotone) curves in 2, we subdivide each curve into a finite number of weakly x-monotone curves, and obtain the set C. Then, we split the set into two disjoint subsets C1 and C2, and we compute their envelope diagrams recursively. Finally, we merge the diagrams in linear time by traversing both diagrams in parallel.

31.2   The Envelope Diagram

The package basically contains two sets of free functions: lower_envelope_x_monotone_2 (begin, end, diag) (similarly upper_envelope_x_monotone_2()) construct the envelope diagram for a given range of x-monotone curves, while lower_envelope_2 (begin, end, diag) (similarly upper_envelope_2()) construct the envelope diagram for a range of arbitrary (not necessarily x-monotone) curves. In this section we explain more on the structure of the envelope diagram these functions output.

The minimization diagram
Figure 31.2:  The lower envelope of eight line segments, labeled A, , H , as constructed in envelope_segments.cpp. The minimization diagram is shown at the bottom, where each diagram vertex points to the point associated with it, and the labels of the segment that induce a diagram edge are displayed below this edge. Note that there exists one edge that represents an overlap (i.e., there are two segments that induce it), and there are also a few edges that represent empty intervals.

A minimization diagram or a maximization diagram is represented by a model of the concept EnvelopeDiagram_1. This concept defines the structure of the subdivision of the x-axis into 0-dimensional cells called vertices, and 1-dimensional cells called edges. The important property of this subdivision is that the identity of the curves that induce the lower envelope (or the upper envelope) over each cell is fixed.

Figure 31.2 shows the lower envelope of a set of eight line segments, and sketches the structure of their minimization diagram. Each diagram vertex v is associated with a point pv on the envelope, which corresponds to either a curve endpoint or to an intersection point of two (or more) curves. The vertex is therefore associated with a set of x-monotone curves that induce the envelope over pv. Each vertex is incident to two edges, one lying to its left and the other to its right.

An edge in the envelope diagram represents a continuous portion of the x-axis, and is associated with a set of x-monotone curves that induce the envelope over this interval. Note that this set may be empty if no x-monotone curves are defined over this interval. In degenerate situations where curves overlap, there may be more than a single curve that induces the envelope over the interval the edge represents. An envelop diagram of a set of curves either consists of a single unbounded edge (in case the curve set is empty or if the envelope contains a single unbounded curve that is below or above all other curves), or at least one vertex and two unbounded edges, while each additional vertex comes along with an additional edge. It is possible to directly access the leftmost edge, representing the unbounded interval that starts at -∞, and the rightmost edge, representing the unbounded interval that ends at . (In the example depicted in Figure 31.2 we have only bounded curves, so the leftmost and rightmost edges represent empty intervals. This is not the case when we deal, for example, with envelopes of sets of lines.)

Any model of the concept must define a geometric traits class, which in turn defines the Point_2 and X_monotone_curve_2 types defined with the diagram features. The geometric traits class must be a model of the ArrangementXMonotoneTraits_2 concept in case we construct envelopes of x-monotone curves. If we are interested in handling arbitrary (not necessarily x-monotone) curves, the traits class must be a model of the ArrangementTraits_2 concept. This concepts refined the ArrangementXMonotoneTraits_2 concept; a traits class that models this concepts must also defines a Curve_2 type, representing an arbitrary planar curve, and provide a functor for subdividing such curves into x-monotone subcurves.

31.3   Examples

31.3.1   Example for Envelope of Line Segments

The following example demonstrates how to compute and traverse the minimization diagram of line segments, as illustrated in Figure 31.2. We use the curve-data traits instantiated by the Arr_segment_traits_2 class, in order to attach a label (a char in this case) to each input segment. We use these labels when we print the minimization diagram:

File: examples/Envelope_2/envelope_segments.cpp
#include <CGAL/Gmpq.h>
#include <CGAL/Cartesian.h>
#include <CGAL/Arr_segment_traits_2.h>
#include <CGAL/Arr_curve_data_traits_2.h>
#include <CGAL/Env_default_diagram_1.h>
#include <CGAL/envelope_2.h>

#include <list>
#include <iostream>

typedef CGAL::Gmpq                                      Number_type;
typedef CGAL::Cartesian<Number_type>                    Kernel;
typedef CGAL::Arr_segment_traits_2<Kernel>              Segment_traits_2;
typedef Segment_traits_2::X_monotone_curve_2            Segment_2;
typedef CGAL::Arr_curve_data_traits_2<Segment_traits_2,
                                      char>             Traits_2;
typedef Traits_2::Point_2                               Point_2;
typedef Traits_2::X_monotone_curve_2                    Labeled_segment_2;
typedef CGAL::Env_default_diagram_1<Traits_2>           Diagram_1;

int main ()
{
  // Consrtuct the input segments and label them 'A' ... 'H'.
  std::list<Labeled_segment_2>   segments;

  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (0, 1),
                                                    Point_2 (2, 3)), 'A'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (1, 2),
                                                    Point_2 (4, 5)), 'B'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (1, 5),
                                                    Point_2 (7, 2)), 'C'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (4, 2),
                                                    Point_2 (6, 4)), 'D'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (8, 3),
                                                    Point_2 (8, 6)), 'E'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (9, 2),
                                                    Point_2 (12, 4)), 'F'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (10, 2),
                                                    Point_2 (12, 1)), 'G'));
  segments.push_back (Labeled_segment_2 (Segment_2 (Point_2 (11, 0),
                                                    Point_2 (11, 5)), 'H'));

  // Compute the minimization diagram that represents their lower envelope.
  Diagram_1              min_diag;

  lower_envelope_x_monotone_2 (segments.begin(), segments.end(),
                               min_diag);

  // Print the minimization diagram.
  Diagram_1::Edge_const_handle     e = min_diag.leftmost();
  Diagram_1::Vertex_const_handle   v;
  Diagram_1::Curve_const_iterator  cit;

  while (e != min_diag.rightmost())
  {
    std::cout << "Edge:";
    if (! e->is_empty())
    {
      for (cit = e->curves_begin(); cit != e->curves_end(); ++cit)
        std::cout << ' ' << cit->data();
    }
    else
      std::cout << " [empty]";
    std::cout << std::endl;

    v = e->right();
    std::cout << "Vertex (" << v->point() << "):";
    for (cit = v->curves_begin(); cit != v->curves_end(); ++cit)
      std::cout << ' ' << cit->data();
    std::cout << std::endl;

    e = v->right();
  }
  CGAL_assertion (e->is_empty());
  std::cout << "Edge: [empty]" << std::endl;

  return (0);
}

31.3.2   Example for Computing the Convex Hull with Envelopes

The next example computes the convex hull of a set of input points by constructing envelopes of unbounded curves, in our case lines that are dual to the input points. Here use the Arr_linear_traits_2 class to compute the lower envelope of the set of dual lines. We read a set of points P = p1, , pn from an input file, and construct the corresponding dual lines P* = p*1, , p*n, where the line p* dual to a point p = (px, py) is given by y = px x - py. We then compute the convex hull of the point-set P, using the fact that the lines that form the lower envelope of P* are dual to the points along the upper part of P's convex hull, and the lines that form the upper envelope of P* are dual to the points along the lower part of the convex hull; see, e.g., [dBvKOS00, Section 11.4] for more details. Note that the leftmost edge of the minimization diagram is associated with the same line as the rightmost edge of the maximization diagram, and vice-verse. We can therefore skip the rightmost edges of both diagrams:

File: examples/Envelope_2/convex_hull.cpp
#include <CGAL/Gmpq.h>
#include <CGAL/Cartesian.h>
#include <CGAL/Arr_linear_traits_2.h>
#include <CGAL/Arr_curve_data_traits_2.h>
#include <CGAL/Env_default_diagram_1.h>
#include <CGAL/envelope_2.h>
#include <vector>

typedef CGAL::Gmpq                                       Number_type;
typedef CGAL::Cartesian<Number_type>                     Kernel;
typedef CGAL::Arr_linear_traits_2<Kernel>                Linear_traits_2;
typedef Linear_traits_2::Point_2                         Point_2;
typedef Linear_traits_2::Line_2                          Line_2;
typedef CGAL::Arr_curve_data_traits_2<Linear_traits_2,
                                      unsigned int>      Traits_2;
typedef Traits_2::X_monotone_curve_2                     Dual_line_2;
typedef CGAL::Env_default_diagram_1<Traits_2>            Diagram_1;

int main (int argc, char **argv)
{
  // Get the name of the input file from the command line, or use the default
  // points.dat file if no command-line parameters are given.
  const char   *filename = "ch_points.dat";

  if (argc > 1)
    filename = argv[1];

  // Open the input file.
  std::ifstream     in_file (filename);

  if (! in_file.is_open())
  {
    std::cerr << "Failed to open " << filename << " ..." << std::endl;
    return (1);
  }

  // Read the points from the file, and construct their dual lines.
  std::list<Dual_line_2>  dual_lines;

  unsigned int            n;
  in_file >> n;
  std::vector<Point_2>    points;
  points.resize (n);

  unsigned int            k;
  for (k = 0; k < n; ++k)
  {
    int px, py;
    in_file >> px >> py;
    points[k] = Point_2 (px, py);

    // The line dual to the point (p_x, p_y) is y = p_x*x - p_y,
    // or: p_x*x - y - p_y = 0:
    Line_2 line = Line_2 (Number_type(px), Number_type(-1), Number_type(-py));

    // Generate the x-monotone curve based on the line and the point index.
    dual_lines.push_back (Dual_line_2 (line, k));
  }
  in_file.close();

  // Compute the lower envelope of dual lines, which corresponds to the upper
  // part of the convex hull, and their upper envelope, which corresponds to
  // the lower part of the convex hull.
  Diagram_1              min_diag;
  Diagram_1              max_diag;

  lower_envelope_x_monotone_2 (dual_lines.begin(), dual_lines.end(),
                               min_diag);

  upper_envelope_x_monotone_2 (dual_lines.begin(), dual_lines.end(),
                               max_diag);

  // Output the points along the boundary convex hull in counterclockwise
  // order. We start by traversing the minimization diagram from left to
  // right, then the maximization diagram from right to left.
  Diagram_1::Edge_const_handle  e = min_diag.leftmost();

  std::cout << "The convex hull of " << points.size() << " input points:";
  while (e != min_diag.rightmost())
  {
    k = e->curve().data();    // The index of the dual point.
    std::cout << " (" << points[k] << ')';
    e = e->right()->right();
  }

  e = max_diag.rightmost();
  while (e != max_diag.leftmost())
  {
    k = e->curve().data();    // The index of the dual point.
    std::cout << " (" << points[k] << ')';
    e = e->left()->left();
  }
  std::cout << std::endl;

  return (0);
}

Envelopes of circles
Figure 31.3:  A set of four circles, as constructed in ex_envelope_circles.cpp. The lower envelope and the upper envelope are shown using thick dashed lines of different colors respectively.

31.3.3   Example for Envelope of Non-Linear Curves

We conclude by an example of envelopes of non-linear curves. We use the Arr_circle_segment_traits_2 class to construct the lower and the upper envelopes of a set of four circles, as depicted in Figure 31.3. Note that unlike the two previous examples, here our curves are not x-monotone, so we use the functions that compute envelopes of arbitrary curves:

File: examples/Envelope_2/envelope_circles.cpp
#include <CGAL/Gmpq.h>
#include <CGAL/Cartesian.h>
#include <CGAL/Arr_circle_segment_traits_2.h>
#include <CGAL/Arrangement_2.h>
#include <CGAL/Env_default_diagram_1.h>
#include <CGAL/envelope_2.h>

typedef CGAL::Gmpq                                    Number_type;
typedef CGAL::Cartesian<Number_type>                  Kernel;
typedef Kernel::Point_2                               Kernel_point_2;
typedef Kernel::Circle_2                              Circle_2;
typedef CGAL::Arr_circle_segment_traits_2<Kernel>     Traits_2;
typedef Traits_2::Curve_2                             Curve_2;
typedef CGAL::Env_default_diagram_1<Traits_2>         Diagram_1;

/*! Print the given envelope diagram. */
void print_diagram (const Diagram_1& diag)
{
  Diagram_1::Edge_const_handle     e = diag.leftmost();
  Diagram_1::Vertex_const_handle   v;

  while (e != diag.rightmost())
  {
    std::cout << "Edge: ";
    if (! e->is_empty())
    {
      Circle_2      circ = e->curve().supporting_circle();
      std::cout << " (x - " << CGAL::to_double(circ.center().x()) << ")^2 +"
                << " (y - " << CGAL::to_double(circ.center().y()) << ")^2 = "
                << CGAL::to_double(circ.squared_radius()) << std::endl;
    }
    else
      std::cout << " [empty]" << std::endl;

    v = e->right();
    std::cout << "Vertex (" << CGAL::to_double(v->point().x()) << ' '
              << CGAL::to_double(v->point().y()) << ')' << std::endl;

    e = v->right();
  }
  CGAL_assertion (e->is_empty());
  std::cout << "Edge: [empty]" << std::endl;

  return;
}

/*! The main program. */
int main ()
{
  // Create four input circles.
  Curve_2          circles[4];

  circles[0] = Circle_2 (Kernel_point_2 (1, 3), CGAL::square(2));
  circles[1] = Circle_2 (Kernel_point_2 (4, 5), CGAL::square(4));
  circles[2] = Circle_2 (Kernel_point_2 (5, 1), CGAL::square(1));
  circles[3] = Circle_2 (Kernel_point_2 (6, 7), CGAL::square(2));

  // Compute the minimization diagram that represents their lower envelope.
  Diagram_1              min_diag;

  lower_envelope_2 (&(circles[0]), &(circles[4]),
                    min_diag);
  print_diagram (min_diag);

  // Compute the maximization diagram that represents the upper envelope.
  Diagram_1              max_diag;

  upper_envelope_2 (&(circles[0]), &(circles[4]),
                    max_diag);
  print_diagram (max_diag);

  return (0);
}