CGAL 4.8.1 - Bounding Volumes
|
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.
CGAL::Min_circle_2<Traits>
CGAL::Min_circle_2_traits_2<K>
MinCircle2Traits
CGAL::Min_ellipse_2<Traits>
CGAL::Min_ellipse_2_traits_2<K>
MinEllipse2Traits
CGAL::Approximate_min_ellipsoid_d<Traits>
ApproximateMinEllipsoid_d_Traits_d
CGAL::min_rectangle_2
CGAL::min_parallelogram_2
CGAL::min_strip_2
CGAL::Min_quadrilateral_default_traits_2<K>
MinQuadrilateralTraits_2
CGAL::rectangular_p_center_2
CGAL::Rectangular_p_center_default_traits_2<K>
RectangularPCenterTraits_2
CGAL::Min_sphere_d<Traits>
CGAL::Min_annulus_d<Traits>
CGAL::Min_sphere_annulus_d_traits_2<K,ET,NT>
CGAL::Min_sphere_annulus_d_traits_3<K,ET,NT>
CGAL::Min_sphere_annulus_d_traits_d<K,ET,NT>
MinSphereAnnulusDTraits
CGAL::Min_sphere_of_spheres_d<Traits>
MinSphereOfSpheresTraits
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... | |
class | 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... | |
class | 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... | |
class | 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... | |
class | 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) |
The function computes a minimum area enclosing parallelogram \( A(P)\) of a given convex point set \( P\). 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) |
The function computes a minimum area enclosing rectangle \( R(P)\) of a given convex point set \( P\). 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) |
The function computes a minimum width enclosing strip \( S(P)\) of a given convex point set \( P\). 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... | |
OutputIterator CGAL::min_parallelogram_2 | ( | ForwardIterator | points_begin, |
ForwardIterator | points_end, | ||
OutputIterator | o, | ||
Traits & | t = Default_traits |
||
) |
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.
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. 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.
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.
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
. OutputIterator
must accept VT
as value type. 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/min_quadrilateral_2.h>
OutputIterator CGAL::min_rectangle_2 | ( | ForwardIterator | points_begin, |
ForwardIterator | points_end, | ||
OutputIterator | o, | ||
Traits & | t = Default_traits |
||
) |
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.
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.
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.
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
. OutputIterator
must accept VT
as value type. 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/min_quadrilateral_2.h>
OutputIterator CGAL::min_strip_2 | ( | ForwardIterator | points_begin, |
ForwardIterator | points_end, | ||
OutputIterator | o, | ||
Traits & | t = Default_traits |
||
) |
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.
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.
If the input range is empty or consists of one element only, o
remains unchanged.
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.
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
. OutputIterator
must accept Traits::Line_2
as value type. 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/min_quadrilateral_2.h>
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 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.
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.
ForwardIterator
must be CGAL::Point_2<K>
for some representation class K
and FT
must be equivalent to K::FT
, Traits
must be a model for RectangularPCenterTraits_2
. OutputIterator
must accept the value type of ForwardIterator
as value type. 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/rectangular_p_center_2.h>