CGAL 5.2 - Bounding Volumes
Bounding Volumes Reference

minCircle.png
Kaspar Fischer, Bernd Gärtner, Thomas Herrmann, Michael Hoffmann, and Sven Schönherr
This package provides algorithms for computing optimal bounding volumes of point sets. In d-dimensional space, the smallest enclosing sphere, ellipsoid (approximate), and annulus can be computed. In 3-dimensional space, the smallest enclosing strip is available as well, and in 2-dimensional space, there are algorithms for a number of additional volumes (rectangles, parallelograms, \( k=2,3,4\) axis-aligned rectangles). The smallest enclosing sphere algorithm can also be applied to a set of d-dimensional spheres.
Introduced in: CGAL 1.1
Depends on: Eigen
BibTeX: cgal:fghhs-bv-20b
License: GPL
Windows Demo: 2D Bounding Volumes
Common Demo Dlls: dlls

Assertions

The optimization code uses infix OPTIMISATION in the assertions, e.g. defining the compiler flag CGAL_OPTIMISATION_NO_PRECONDITIONS switches precondition checking off, cf. Section Checks.

Classified Reference Pages

Bounding Areas and Volumes

Modules

 Concepts
 

Classes

class  CGAL::Approximate_min_ellipsoid_d< Traits >
 An object of class Approximate_min_ellipsoid_d is an approximation to the ellipsoid of smallest volume enclosing a finite multiset of points in \( d\)-dimensional Euclidean space \( \E^d\), \( d\ge 2\). More...
 
struct  CGAL::Approximate_min_ellipsoid_d_traits_2< K, ET >
 The class Approximate_min_ellipsoid_d_traits_2 is a traits class for CGAL::Approximate_min_ellipsoid_d<Traits> using the 2-dimensional CGAL kernel. More...
 
struct  CGAL::Approximate_min_ellipsoid_d_traits_3< K, ET >
 The class Approximate_min_ellipsoid_d_traits_3 is a traits class for CGAL::Approximate_min_ellipsoid_d<Traits> using the 3-dimensional CGAL kernel. More...
 
struct  CGAL::Approximate_min_ellipsoid_d_traits_d< K, ET >
 The class Approximate_min_ellipsoid_d_traits_d is a traits class for CGAL::Approximate_min_ellipsoid_d<Traits> using the d-dimensional CGAL kernel. More...
 
class  CGAL::Min_annulus_d< Traits >
 An object of the class Min_annulus_d is the unique annulus (region between two concentric spheres with radii \( r\) and \( R\), \( r \leq R\)) enclosing a finite set of points in \( d\)-dimensional Euclidean space \( \E^d\), where the difference \( R^2-r^2\) is minimal. More...
 
class  CGAL::Min_circle_2< Traits >
 An object of the class Min_circle_2 is the unique circle of smallest area enclosing a finite (multi)set of points in two-dimensional Euclidean space \( \E^2\). More...
 
class  CGAL::Min_circle_2_traits_2< K >
 The class Min_circle_2_traits_2 is a traits class for Min_circle_2<Traits> using the two-dimensional CGAL kernel. More...
 
class  CGAL::Min_ellipse_2< Traits >
 An object of the class Min_ellipse_2 is the unique ellipse of smallest area enclosing a finite (multi)set of points in two-dimensional euclidean space \( \E^2\). More...
 
class  CGAL::Min_ellipse_2_traits_2< K >
 The class Min_ellipse_2_traits_2 is a traits class for CGAL::Min_ellipse_2<Traits> using the two-di-men-sional CGAL kernel. More...
 
struct  CGAL::Min_quadrilateral_default_traits_2< K >
 The class Min_quadrilateral_default_traits_2 is a traits class for the functions min_rectangle_2, min_parallelogram_2 and min_strip_2 using a two-dimensional CGAL kernel. More...
 
class  CGAL::Min_sphere_annulus_d_traits_2< K, ET, NT >
 The class Min_sphere_annulus_d_traits_2 is a traits class for the \( d\)-dimensional optimisation algorithms using the two-dimensional CGAL kernel. More...
 
class  CGAL::Min_sphere_annulus_d_traits_3< K, ET, NT >
 The class Min_sphere_annulus_d_traits_3 is a traits class for the \( d\)-dimensional optimisation algorithms using the three-dimensional CGAL kernel. More...
 
class  CGAL::Min_sphere_annulus_d_traits_d< K, ET, NT >
 The class Min_sphere_annulus_d_traits_d is a traits class for the \( d\)-dimensional optimisation algorithms using the \( d\)-dimensional CGAL kernel. More...
 
class  CGAL::Min_sphere_d< Traits >
 An object of the class Min_sphere_d is the unique sphere of smallest volume enclosing a finite (multi)set of points in \( d\)-dimensional Euclidean space \( \E^d\). More...
 
class  CGAL::Min_sphere_of_points_d_traits_2< K, FT, UseSqrt, Algorithm >
 The class Min_sphere_of_points_d_traits_2<K,FT,UseSqrt,Algorithm> is a model for concept MinSphereOfSpheresTraits. More...
 
class  CGAL::Min_sphere_of_points_d_traits_3< K, FT, UseSqrt, Algorithm >
 The class Min_sphere_of_points_d_traits_3<K,FT,UseSqrt,Algorithm> is a model for concept MinSphereOfSpheresTraits. More...
 
class  CGAL::Min_sphere_of_points_d_traits_d< K, FT, Dim, UseSqrt, Algorithm >
 The class Min_sphere_of_points_d_traits_d<K,FT,Dim,UseSqrt,Algorithm> is a model for concept MinSphereOfSpheresTraits. More...
 
class  CGAL::Min_sphere_of_spheres_d< Traits >
 An object of the class Min_sphere_of_spheres_d is a data structure that represents the unique sphere of smallest volume enclosing a finite set of spheres in \( d\)-dimensional Euclidean space \( \E^d\). More...
 
class  CGAL::Min_sphere_of_spheres_d_traits_2< K, FT, UseSqrt, Algorithm >
 The class Min_sphere_of_spheres_d_traits_2<K,FT,UseSqrt,Algorithm> is a model for concept MinSphereOfSpheresTraits. More...
 
class  CGAL::Min_sphere_of_spheres_d_traits_3< K, FT, UseSqrt, Algorithm >
 The class Min_sphere_of_spheres_d_traits_3<K,FT,UseSqrt,Algorithm> is a model for concept MinSphereOfSpheresTraits. More...
 
class  CGAL::Min_sphere_of_spheres_d_traits_d< K, FT, Dim, UseSqrt, Algorithm >
 The class Min_sphere_of_spheres_d_traits_d<K,FT,Dim,UseSqrt,Algorithm> is a model for concept MinSphereOfSpheresTraits. More...
 
class  CGAL::Rectangular_p_center_default_traits_2< K >
 The class Rectangular_p_center_default_traits_2 defines types and operations needed to compute rectilinear \( p\)-centers of a planar point set using the function rectangular_p_center_2(). More...
 

Functions

template<class ForwardIterator , class OutputIterator , class Traits >
OutputIterator CGAL::min_parallelogram_2 (ForwardIterator points_begin, ForwardIterator points_end, OutputIterator o, Traits &t=Default_traits)
 computes a minimum area enclosing parallelogram 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. More...
 
template<class ForwardIterator , class OutputIterator , class Traits >
OutputIterator CGAL::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. More...
 
template<class ForwardIterator , class OutputIterator , class Traits >
OutputIterator CGAL::min_strip_2 (ForwardIterator points_begin, ForwardIterator points_end, OutputIterator o, Traits &t=Default_traits)
 computes a minimum enclosing strip of the point set described by [points_begin, points_end), writes its two bounding lines to o and returns the past-the-end iterator of this sequence. More...
 
template<class ForwardIterator , class OutputIterator , class FT , class Traits >
OutputIterator CGAL::rectangular_p_center_2 (ForwardIterator f, ForwardIterator l, OutputIterator o, FT &r, int p, const Traits &t=Default_traits)
 Computes rectilinear \( p\)-centers of a planar point set, i.e. a set of \( p\) points such that the maximum minimal \( L_{\infty}\)-distance between both sets is minimized. More...
 

Function Documentation

◆ min_parallelogram_2()

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

#include <CGAL/min_quadrilateral_2.h>

computes a minimum area enclosing parallelogram 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.

The function computes a minimum area enclosing parallelogram \( A(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 parallelogram 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.

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.

  1. If Traits is specified, it must be a model for MinQuadrilateralTraits_2 and the value type VT of ForwardIterator is Traits::Point_2. Otherwise VT must be CGAL::Point_2<K> for some kernel K.
  2. OutputIterator must accept VT as value type.
See also
CGAL::min_rectangle_2()
CGAL::min_strip_2()
MinQuadrilateralTraits_2
CGAL::Min_quadrilateral_default_traits_2<K>

Implementation

We use a rotating caliper algorithm [12], [15] 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 parallelogram of P.


File Min_quadrilateral_2/minimum_enclosing_parallelogram_2.cpp

#include <CGAL/Simple_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>
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 parallelogram p_m of p
Polygon_2 p_m;
p.vertices_begin(), p.vertices_end(), std::back_inserter(p_m));
std::cout << p_m << std::endl;
return 0;
}
Examples:
Min_quadrilateral_2/minimum_enclosing_parallelogram_2.cpp.

◆ min_rectangle_2()

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

#include <CGAL/min_quadrilateral_2.h>

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.

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.

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.

  1. If Traits is specified, it must be a model for MinQuadrilateralTraits_2 and the value type VT of ForwardIterator is Traits::Point_2. Otherwise VT must be CGAL::Point_2<K> for some kernel K.
  2. OutputIterator must accept 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 [14] 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 Min_quadrilateral_2/minimum_enclosing_rectangle_2.cpp

#include <CGAL/Simple_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>
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;
p.vertices_begin(), p.vertices_end(), std::back_inserter(p_m));
std::cout << p_m << std::endl;
return 0;
}
Examples:
Min_quadrilateral_2/minimum_enclosing_rectangle_2.cpp.

◆ min_strip_2()

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

#include <CGAL/min_quadrilateral_2.h>

computes a minimum enclosing strip of the point set described by [points_begin, points_end), writes its two bounding lines to o and returns the past-the-end iterator of this sequence.

The function computes a minimum width enclosing strip \( S(P)\) of a given convex point set \( P\). A strip is the closed region bounded by two parallel lines in the plane. Note that \( S(P)\) is not unique in general. The focus on convex sets is no restriction, since any parallelogram 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.

If the input range is empty or consists of one element only, o remains unchanged.

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.

  1. If Traits is specified, it must be a model for MinQuadrilateralTraits_2 and the value type VT of ForwardIterator is Traits::Point_2. Otherwise VT must be CGAL::Point_2<K> for some kernel K.
  2. OutputIterator must accept Traits::Line_2 as value type.
See also
CGAL::min_rectangle_2()
CGAL::min_parallelogram_2()
MinQuadrilateralTraits_2
CGAL::Min_quadrilateral_default_traits_2<K>

Implementation

We use a rotating caliper algorithm [14] 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 strip of P.


File Min_quadrilateral_2/minimum_enclosing_strip_2.cpp

#include <CGAL/Simple_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>
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 strip p_m of p
Line_2 p_m[2];
CGAL::min_strip_2(p.vertices_begin(), p.vertices_end(), p_m);
std::cout << p_m[0] << "\n" << p_m[1] << std::endl;
return 0;
}
Examples:
Min_quadrilateral_2/minimum_enclosing_strip_2.cpp.

◆ rectangular_p_center_2()

template<class ForwardIterator , class OutputIterator , class FT , class Traits >
OutputIterator CGAL::rectangular_p_center_2 ( ForwardIterator  f,
ForwardIterator  l,
OutputIterator  o,
FT &  r,
int  p,
const Traits &  t = Default_traits 
)

#include <CGAL/rectangular_p_center_2.h>

Computes rectilinear \( p\)-centers of a planar point set, i.e. a set of \( p\) points such that the maximum minimal \( L_{\infty}\)-distance between both sets is minimized.

More formally the problem can be defined as follows.

Given a finite set \( \mathcal{P}\) of points, compute a point set \( \mathcal{C}\) with \( |\mathcal{C}| \le p\) such that the \( p\)-radius of \( \mathcal{P}\),

\[ rad_p(\mathcal{P}) := \max_{P \in \mathcal{P}} \min_{Q \in \mathcal{C}} || P - Q ||_\infty \]

is minimized. We can interpret \( \mathcal{C}\) as the best approximation (with respect to the given metric) for \( \mathcal{P}\) with at most \( p\) points.

computes rectilinear p-centers for the point set described by the range [f, l), sets r to the corresponding \( p\)-radius, writes the at most p center points to o and returns the past-the-end iterator of this sequence.

Precondition
2 \( \le\) p \( \le\) 4.

The geometric types and operations to be used for the computation are specified by the traits class parameter t. This parameter can be omitted if ForwardIterator refers to a point type from the 2D-Kernel. In this case, a default traits class (Rectangular_p_center_default_traits_2<K>) is used.

  1. Either: (if no traits parameter is given) Value type of ForwardIterator must be CGAL::Point_2<K> for some representation class K and FT must be equivalent to K::FT,
  2. Or: (if a traits parameter is specified) Traits must be a model for RectangularPCenterTraits_2.
  3. OutputIterator must accept the value type of ForwardIterator as value type.
See also
RectangularPCenterTraits_2
CGAL::Rectangular_p_center_default_traits_2<K>
CGAL::sorted_matrix_search()

Implementation

The runtime is linear for \( p \in \{2,\,3\}\) and \( \mathcal{O}(n \cdot \log n)\) for \( p = 4\) where \( n\) is the number of input points. These runtimes are worst case optimal. The \( 3\)-center algorithm uses a prune-and-search technique described in [9]. The \( 4\)-center implementation uses sorted matrix search [1], [2] and fast algorithms for piercing rectangles [13].

Example

The following code generates a random set of ten points and computes its two-centers.


File Rectangular_p_center_2/rectangular_p_center_2.cpp

#include <CGAL/Simple_cartesian.h>
#include <CGAL/point_generators_2.h>
#include <CGAL/rectangular_p_center_2.h>
#include <CGAL/IO/Ostream_iterator.h>
#include <CGAL/algorithm.h>
#include <iostream>
#include <algorithm>
#include <vector>
typedef double FT;
typedef Kernel::Point_2 Point;
typedef std::vector<Point> Cont;
typedef CGAL::Random_points_in_square_2<Point> Generator;
int main()
{
int n = 10;
int p = 2;
OIterator cout_ip(std::cout);
Cont points;
std::copy_n(Generator(1), n, std::back_inserter(points));
std::cout << "Generated Point Set:\n";
std::copy(points.begin(), points.end(), cout_ip);
FT p_radius;
std::cout << "\n\n" << p << "-centers:\n";
points.begin(), points.end(), cout_ip, p_radius, 3);
std::cout << "\n\n" << p << "-radius = " << p_radius << std::endl;
return 0;
}
Examples:
Rectangular_p_center_2/rectangular_p_center_2.cpp.