CGAL::min_rectangle_2

Definition

The function computes a minimum area enclosing rectangle R(P) of a given convex point set P. Note that R(P) is not necessarily axis-parallel, and it is in general not unique. The focus on convex sets is no restriction, since any rectangle enclosing P - as a convex set - contains the convex hull of P. For general point sets one has to compute the convex hull as a preprocessing step.

#include <CGAL/min_quadrilateral_2.h>

template < class ForwardIterator, class OutputIterator, class Traits >
OutputIterator
min_rectangle_2 (
ForwardIterator points_begin,
ForwardIterator points_end,
OutputIterator o,
Traits& t = Default_traits)

computes a minimum area enclosing rectangle of the point set described by [points_begin, points_end), writes its vertices (counterclockwise) to o, and returns the past-the-end iterator of this sequence.
If the input range is empty, o remains unchanged.
If the input range consists of one element only, this point is written to o four times.

Precondition

The points denoted by the range [points_begin, points_end) form the boundary of a simple convex polygon P in counterclockwise orientation.

The geometric types and operations to be used for the computation are specified by the traits class parameter t. The parameter can be omitted, if ForwardIterator refers to a two-dimensional point type from one the Cgal kernels. In this case, a default traits class (Min_quadrilateral_default_traits_2<K>) is used.

Requirement

  1. If Traits is specified, it is a model for MinQuadrilateralTraits_2 and the value type VT of ForwardIterator is Traits::Point_2. Otherwise VT is CGAL::Point_2<K> for some kernel K.
  2. OutputIterator accepts VT as value type.

See Also

CGAL::min_parallelogram_2
CGAL::min_strip_2
MinQuadrilateralTraits_2
CGAL::Min_quadrilateral_default_traits_2<K>

Implementation

We use a rotating caliper algorithm [Tou83] with worst case running time linear in the number of input points.

Example

The following code generates a random convex polygon P with 20 vertices and computes the minimum enclosing rectangle of P.

File: examples/Min_quadrilateral_2/minimum_enclosing_rectangle_2.cpp
#include <CGAL/Cartesian.h>
#include <CGAL/Polygon_2.h>
#include <CGAL/point_generators_2.h>
#include <CGAL/random_convex_set_2.h>
#include <CGAL/min_quadrilateral_2.h>
#include <iostream>

struct Kernel : public CGAL::Cartesian<double> {};

typedef Kernel::Point_2                           Point_2;
typedef Kernel::Line_2                            Line_2;
typedef CGAL::Polygon_2<Kernel>                   Polygon_2;
typedef CGAL::Random_points_in_square_2<Point_2>  Generator;

int main()
{
  // build a random convex 20-gon p
  Polygon_2 p;
  CGAL::random_convex_set_2(20, std::back_inserter(p), Generator(1.0));
  std::cout << p << std::endl;

  // compute the minimal enclosing rectangle p_m of p
  Polygon_2 p_m;
  CGAL::min_rectangle_2(
    p.vertices_begin(), p.vertices_end(), std::back_inserter(p_m));
  std::cout << p_m << std::endl;

  return 0;
}