An object of the class Min_annulus_d<Traits> is the unique annulus (region between two concentric spheres with radii and , ) enclosing a finite set of points in -dimensional Euclidean space , where the difference is minimal. For a point set we denote by the smallest annulus that contains all points of . Note that can be degenerate, i.e. Ø if Ø and if .
An inclusion-minimal subset of with is called a support set, the points in are the support points. A support set has size at most , and all its points lie on the boundary of . In general, the support set is not necessarily unique.
The underlying algorithm can cope with all kinds of input, e.g. may be empty or points may occur more than once. The algorithm computes a support set which remains fixed until the next set, insert, or clear operation.
#include <CGAL/Min_annulus_d.h>
The template parameter Traits is a model for OptimisationDTraits.
We provide the models Optimisation_d_traits_2, Optimisation_d_traits_3, and Optimisation_d_traits_d using the two-, three-, and -dimensional CGAL kernel, respectively.
| |
typedef to Traits::Point_d.
Point type used to represent the input points.
| |
| |
typedef to Traits::FT.
Number type used to return the squared radii of the smallest
enclosing annulus.
| |
| |
typedef to Traits::ET.
Number type used to do the exact computations in the underlying
solver for quadratic programs (cf. Implementation).
| |
| |
non-mutable model of the STL concept RandomAccessIterator
with value type Point. Used to access the points
of the smallest enclosing annulus.
| |
| |
non-mutable model of the STL concept RandomAccessIterator
with value type Point. Used to access the support points
of the smallest enclosing annulus.
| |
| |
non-mutable model of the STL concept RandomAccessIterator
with value type Point. Used to access the inner support points
of the smallest enclosing annulus.
| |
| |
non-mutable model of the STL concept RandomAccessIterator
with value type Point. Used to access the outer support points
of the smallest enclosing annulus.
| |
| |
non-mutable model of the STL concept RandomAccessIterator
with value type ET. Used to access the coordinates of
the center of the smallest enclosing annulus.
|
| |||||
initializes min_annulus to Ø.
| |||||
| |||||
| |||||
initializes min_annulus to with being the set of points
in the range [first,last).
|
|
| returns the dimension of the points in . If min_annulus is empty, the ambient dimension is . |
|
| returns the number of points of min_annulus, i.e. . |
|
| |
returns the number of support points of min_annulus, i.e. . | ||
|
| |
returns the number of support points of min_annulus which lie on the inner sphere. | ||
|
| |
returns the number of support points of min_annulus which lie on the outer sphere. |
|
| returns an iterator referring to the first point of min_annulus. |
|
| returns the corresponding past-the-end iterator. |
|
| |
returns an iterator referring to the first inner support point of min_annulus. | ||
|
| |
returns the corresponding past-the-end iterator. |
|
| |||||
returns an iterator referring to the first outer support point of min_annulus. | ||||||
|
| |||||
returns the corresponding past-the-end iterator. | ||||||
|
|
returns the center of min_annulus.
| ||||
|
| |||||
returns the squared inner radius of min_annulus.
| ||||||
|
| |||||
returns the squared outer radius of min_annulus.
|
|
| |||
returns an iterator referring to the first coordinate
of the center of min_annulus.
| ||||
|
| |||
returns the corresponding past-the-end iterator. | ||||
|
| |||
returns the numerator of the squared inner radius of min_annulus. | ||||
|
| |||
returns the numerator of the squared outer radius of min_annulus. | ||||
|
| |||
returns the denominator of the squared radii of min_annulus. |
The bounded area of the smallest enclosing annulus lies between the inner and the outer sphere. The boundary is the union of both spheres. By definition, an empty annulus has no boundary and no bounded side, i.e. its unbounded side equals the whole space .
|
| |||
returns CGAL::ON_BOUNDED_SIDE,
CGAL::ON_BOUNDARY, or
CGAL::ON_UNBOUNDED_SIDE iff p lies
properly inside, on the boundary, or properly outside of
min_annulus, resp.
| ||||
|
| |||
returns true, iff p lies properly inside min_annulus.
| ||||
|
| |||
returns true, iff p lies on the boundary of min_annulus.
| ||||
|
| |||
returns true, iff p lies properly outside of min_annulus.
| ||||
|
| returns true, iff min_annulus is empty (this implies degeneracy). | ||
|
| returns true, iff min_annulus is degenerate, i.e. if min_annulus is empty or equal to a single point. |
|
| resets min_annulus to Ø. | ||||
| ||||||
|
| |||||
sets min_annulus to , where is the set of points in
the range [first,last).
| ||||||
|
|
inserts p into min_annulus.
| ||||
| ||||||
|
| |||||
inserts the points in the range [first,last) into
min_annulus and recomputes the smallest enclosing annulus.
|
An object min_annulus is valid, iff
|
| |
returns true, iff min_annulus is valid. If verbose is true, some messages concerning the performed checks are written to standard error stream. The second parameter level is not used, we provide it only for consistency with interfaces of other classes. |
|
| returns a const reference to the traits class object. |
|
|
writes min_annulus to output stream os.
|
|
|
reads min_annulus from input stream is.
|
CGAL::Min_sphere_d<Traits>
CGAL::Optimisation_d_traits_2<K,ET,NT>
CGAL::Optimisation_d_traits_3<K,ET,NT>
CGAL::Optimisation_d_traits_d<K,ET,NT>
OptimisationDTraits
The problem of finding the smallest enclosing annulus of a finite point set can be formulated as an optimization problem with linear constraints and a linear objective function. The solution is obtained using our exact solver for linear and quadratic programs [GS00].
The creation time is almost always linear in the number of points. Access functions and predicates take constant time, inserting a point takes almost always linear time. The clear operation and the check for validity each take linear time.