CGAL::Min_annulus_d<Traits>

Definition

An object of the class Min_annulus_d<Traits> is the unique annulus (region between two concentric spheres with radii r and R, r R) enclosing a finite set of points in d-dimensional Euclidean space d, where the difference R2-r2 is minimal. For a point set P we denote by ma(P) the smallest annulus that contains all points of P. Note that ma(P) can be degenerate, i.e. ma(P)=Ø if P=Ø and ma(P)={p} if P={p}.

An inclusion-minimal subset S of P with ma(S)=ma(P) is called a support set, the points in S are the support points. A support set has size at most d+2, and all its points lie on the boundary of ma(P). In general, the support set is not necessarily unique.

The underlying algorithm can cope with all kinds of input, e.g. P may be empty or points may occur more than once. The algorithm computes a support set S which remains fixed until the next set, insert, or clear operation.

#include <CGAL/Min_annulus_d.h>

Requirements

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 d-dimensional CGAL kernel, respectively.

Types

Min_annulus_d<Traits>::Point
typedef to Traits::Point_d. Point type used to represent the input points.


Min_annulus_d<Traits>::FT
typedef to Traits::FT. Number type used to return the squared radii of the smallest enclosing annulus.


Min_annulus_d<Traits>::ET
typedef to Traits::ET. Number type used to do the exact computations in the underlying solver for quadratic programs (cf. Implementation).


Min_annulus_d<Traits>::Point_iterator
non-mutable model of the STL concept RandomAccessIterator with value type Point. Used to access the points of the smallest enclosing annulus.


Min_annulus_d<Traits>::Support_point_iterator
non-mutable model of the STL concept RandomAccessIterator with value type Point. Used to access the support points of the smallest enclosing annulus.


Min_annulus_d<Traits>::Inner_support_point_iterator
non-mutable model of the STL concept RandomAccessIterator with value type Point. Used to access the inner support points of the smallest enclosing annulus.


Min_annulus_d<Traits>::Outer_support_point_iterator
non-mutable model of the STL concept RandomAccessIterator with value type Point. Used to access the outer support points of the smallest enclosing annulus.


Min_annulus_d<Traits>::Coordinate_iterator
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.

Creation

Min_annulus_d<Traits> min_annulus ( Traits traits = Traits(),
int verbose = 0,
std::ostream& stream = std::cout);
initializes min_annulus to ma(Ø).


template < class InputIterator >
Min_annulus_d<Traits> min_annulus ( InputIterator first,
InputIterator last,
Traits traits = Traits(),
int verbose = 0,
std::ostream& stream = std::cout);
initializes min_annulus to ma(P) with P being the set of points in the range [first,last).
Requirement: The value type of InputIterator is Point.
Precondition: All points have the same dimension.


begin of advanced section  advanced  begin of advanced section
If verbose is set to 1, 2, or 3 then some, more, or full verbose output of the underlying solver for linear programs is written to stream, resp.
end of advanced section  advanced  end of advanced section

Access Functions

int min_annulus.ambient_dimension ()
returns the dimension of the points in P. If min_annulus is empty, the ambient dimension is -1.

int min_annulus.number_of_points ()
returns the number of points of min_annulus, i.e. |P|.

int min_annulus.number_of_support_points ()
returns the number of support points of min_annulus, i.e. |S|.

int min_annulus.number_of_inner_support_points ()
returns the number of support points of min_annulus which lie on the inner sphere.

int min_annulus.number_of_outer_support_points ()
returns the number of support points of min_annulus which lie on the outer sphere.



Point_iterator min_annulus.points_begin ()
returns an iterator referring to the first point of min_annulus.
Point_iterator min_annulus.points_end ()
returns the corresponding past-the-end iterator.



Support_point_iterator
min_annulus.support_points_begin ()
returns an iterator referring to the first support point of min_annulus.
Support_point_iterator
min_annulus.support_points_end ()
returns the corresponding past-the-end iterator.



Inner_support_point_iterator
min_annulus.inner_support_points_begin ()
returns an iterator referring to the first inner support point of min_annulus.
Inner_support_point_iterator
min_annulus.inner_support_points_end ()
returns the corresponding past-the-end iterator.



Outer_support_point_iterator
min_annulus.outer_support_points_begin ()
returns an iterator referring to the first outer support point of min_annulus.
Outer_support_point_iterator
min_annulus.outer_support_points_end ()
returns the corresponding past-the-end iterator.

Point min_annulus.center ()
returns the center of min_annulus.
Requirement: An implicit conversion from ET to RT is available.
Precondition: min_annulus is not empty.

FT min_annulus.squared_inner_radius ()
returns the squared inner radius of min_annulus.
Requirement: An implicit conversion from ET to RT is available.
Precondition: min_annulus is not empty.

FT min_annulus.squared_outer_radius ()
returns the squared outer radius of min_annulus.
Requirement: An implicit conversion from ET to RT is available.
Precondition: min_annulus is not empty.



Coordinate_iterator
min_annulus.center_coordinates_begin ()
returns an iterator referring to the first coordinate of the center of min_annulus.
Note: The coordinates have a rational representation, i.e. the first d elements of the iterator range are the numerators and the (d+1)-st element is the common denominator.
Coordinate_iterator
min_annulus.center_coordinates_end ()
returns the corresponding past-the-end iterator.

ET min_annulus.squared_inner_radius_numerator ()
returns the numerator of the squared inner radius of min_annulus.

ET min_annulus.squared_outer_radius_numerator ()
returns the numerator of the squared outer radius of min_annulus.

ET min_annulus.squared_radii_denominator ()
returns the denominator of the squared radii of min_annulus.

Predicates

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 d.

CGAL::Bounded_side min_annulus.bounded_side ( Point p)
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.
Precondition: The dimension of p equals min_annulus.ambient_dimension() if min_annulus is not empty.

bool min_annulus.has_on_bounded_side ( Point p)
returns true, iff p lies properly inside min_annulus.
Precondition: The dimension of p equals min_annulus.ambient_dimension() if min_annulus is not empty.

bool min_annulus.has_on_boundary ( Point p)
returns true, iff p lies on the boundary of min_annulus.
Precondition: The dimension of p equals min_annulus.ambient_dimension() if min_annulus is not empty.

bool min_annulus.has_on_unbounded_side ( Point p)
returns true, iff p lies properly outside of min_annulus.
Precondition: The dimension of p equals min_annulus.ambient_dimension() if min_annulus is not empty.

bool min_annulus.is_empty ()
returns true, iff min_annulus is empty (this implies degeneracy).

bool min_annulus.is_degenerate ()
returns true, iff min_annulus is degenerate, i.e. if min_annulus is empty or equal to a single point.

Modifiers

void min_annulus.clear ()
resets min_annulus to ma(Ø).

template < class InputIterator >
void
min_annulus.set ( InputIterator first,
InputIterator last)
sets min_annulus to ma(P), where P is the set of points in the range [first,last).
Requirement: The value type of InputIterator is Point.
Precondition: All points have the same dimension.

void min_annulus.insert ( Point p)
inserts p into min_annulus.
Precondition: The dimension of p equals min_annulus.ambient_dimension() if min_annulus is not empty.

template < class InputIterator >
void
min_annulus.insert ( InputIterator first,
InputIterator last)
inserts the points in the range [first,last) into min_annulus and recomputes the smallest enclosing annulus.
Requirement: The value type of InputIterator is Point.
Precondition: All points have the same dimension. If min_annulus is not empty, this dimension must be equal to min_annulus.ambient_dimension().

Validity Check

An object min_annulus is valid, iff

Note: In this release only the first item is considered by the validity check.

bool
min_annulus.is_valid ( bool verbose = false,
int level = 0)
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.

Miscellaneous

const Traits& min_annulus.traits ()
returns a const reference to the traits class object.

I/O

std::ostream& std::ostream& os << min_annulus
writes min_annulus to output stream os.
Requirement: The output operator is defined for Point.

std::istream& std::istream& is >> min_annulus&
reads min_annulus from input stream is.
Requirement: The input operator is defined for Point.

See Also

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

Implementation

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.