Ron Wein
A continuous curve in is called -monotone, if every vertical line intersects it at a single point at most. For example, the circle is not -monotone as the vertical line intersects it at and at ; however, it is possible to split the circle into an upper part and a lower part, such that both of these parts are -monotone. We consider vertical segments as weakly -monotone, to properly handle inputs that contain such vertical curves.
An -monotone curve can be represented as a univariate function , defined over some continuous range . Given a set of -monotone curves, their lower envelope is defined as the point-wise minimum of all curves. Namely, the lower envelope of the set can be defined as the following function:
where we define
Similarly, the upper envelope of
(x) = max
where in this case
Given a set of
Lower and upper envelopes can be efficiently computed using a
divide-and-conquer approach. First, note that the envelope diagram for
a single
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
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
Figure 23.2 shows the lower envelope of a set of
eight line segments, and sketches the structure of their minimization
diagram. Each diagram vertex
An edge in the envelope diagram represents a continuous portion of the
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
The following example demonstrates how to compute and traverse the minimization diagram of line segments, as illustrated in Figure 23.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/ex_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); }
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
File: examples/Envelope_2/ex_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. unsigned int n; int px, py; std::vector<Point_2> points; Line_2 line; std::list<Dual_line_2> dual_lines; unsigned int k; in_file >> n; points.resize (n); for (k = 0; k < n; ++k) { 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 = 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); }
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 23.3. Note that unlike the two previous
examples, here our curves are not
File: examples/Envelope_2/ex_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); }