An object of the class Min_sphere_d<Traits> is the unique sphere of smallest volume enclosing a finite (multi)set of points in d-dimensional Euclidean space d. For a set P we denote by ms(P) the smallest sphere that contains all points of P. ms(P) can be degenerate, i.e. ms(P)=∅ if P=∅ and ms(P)={p} if P={p}.
An inclusion-minimal subset S of P with ms(S)=ms(P) is called a support set, the points in S are the support points. A support set has size at most d+1, and all its points lie on the boundary of ms(P). In general, neither the support set nor its size are unique.
The algorithm computes a support set S which remains fixed until the next insert or clear operation.
Please note: This class is (almost) obsolete. The class CGAL::Min_sphere_of_spheres_d<Traits> solves a more general problem and is faster then Min_sphere_d<Traits> even if used only for points as input. Most importantly, CGAL::Min_sphere_of_spheres_d<Traits> has a specialized implementation for floating-point arithmetic which ensures correct results in a large number of cases (including highly degenerate ones). In contrast, Min_sphere_d<Traits> is not reliable under floating-point computations. The only advantage of Min_sphere_d<Traits> over CGAL::Min_sphere_of_spheres_d<Traits> is that the former can deal with points in homogeneous coordinates, in which case the algorithm is division-free. Thus, Min_sphere_d<Traits> might still be an option in case your input number type cannot (efficiently) divide.
#include <CGAL/Min_sphere_d.h>
The class Min_sphere_d<Traits> expects a model of the concept OptimisationDTraits as its template argument. We provide the models CGAL::Optimisation_d_traits_2, CGAL::Optimisation_d_traits_3 and CGAL::Optimisation_d_traits_d for two-, three-, and d-dimensional points respectively.
| |||||
creates a variable of type Min_sphere_d<Traits> and
initializes it to ms(∅).
If the traits parameter is not supplied, the class Traits
must provide a default constructor.
| |||||
| |||||
| |||||
creates a variable min_sphere of type Min_sphere_d<Traits>.
It is initialized to ms(P) with P being the set of points
in the range [first,last).
|
|
| |||
returns the number of points of min_sphere, i.e. |P|. | ||||
|
| |||
returns the number of support points of min_sphere, i.e. |S|. | ||||
|
| returns an iterator referring to the first point of min_sphere. | ||
|
| returns the corresponding past-the-end iterator. | ||
|
| |||
returns an iterator referring to the first support point of min_sphere. | ||||
|
| |||
returns the corresponding past-the-end iterator. | ||||
|
| |||
returns the dimension of the points in P. If min_sphere is empty, the ambient dimension is -1. | ||||
|
|
returns the center of min_sphere.
| ||
|
| |||
returns the squared radius of min_sphere.
|
By definition, an empty Min_sphere_d<Traits> has no boundary and no bounded side, i.e. its unbounded side equals the whole space d.
|
| |||
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_sphere, resp.
| ||||
|
| |||
returns true, iff p lies properly inside min_sphere.
| ||||
|
| |||
returns true, iff p lies on the boundary
of min_sphere.
| ||||
|
| |||
returns true, iff p lies properly outside of min_sphere.
| ||||
|
| returns true, iff min_sphere is empty (this implies degeneracy). | ||
|
| returns true, iff min_sphere is degenerate, i.e. if min_sphere is empty or equal to a single point, equivalently if the number of support points is less than 2. |
|
| resets min_sphere to ms(∅). | ||||
| ||||||
|
| |||||
sets min_sphere to the ms(P), where P is the set of points
in the range [first,last).
| ||||||
|
|
inserts p into min_sphere. If p lies inside the
current sphere, this is a constant-time operation, otherwise
it might take longer, but usually substantially less than
recomputing the smallest enclosing sphere from scratch.
| ||||
| ||||||
|
| |||||
inserts the points in the range [first,last)
into min_sphere and recomputes the smallest enclosing sphere, by
calling insert for all points in the range.
|
Note: Under inexact arithmetic, the result of the validation is not realiable, because the checker itself can suffer from numerical problems.
|
| |
returns true, iff min_sphere 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_sphere to output stream os.
|
|
|
reads min_sphere from input stream is.
|
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
CGAL::Min_circle_2<Traits>
CGAL::Min_sphere_of_spheres_d<Traits>
CGAL::Min_annulus_d<Traits>
We implement the algorithm of Welzl with move-to-front heuristic [Wel91] for small point sets, combined with a new efficient method for large sets, which is particularly tuned for moderately large dimension (d ≤ 20) [Gär99]. The creation time is almost always linear in the number of points. Access functions and predicates take constant time, inserting a point might take up to linear time, but substantially less than computing the new smallest enclosing sphere from scratch. The clear operation and the check for validity each take linear time.
File: examples/Min_sphere_d/min_sphere_d.cpp
#include <CGAL/Cartesian_d.h> #include <iostream> #include <cstdlib> #include <CGAL/Random.h> #include <CGAL/Min_sphere_annulus_d_traits_d.h> #include <CGAL/Min_sphere_d.h> typedef CGAL::Cartesian_d<double> K; typedef CGAL::Min_sphere_annulus_d_traits_d<K> Traits; typedef CGAL::Min_sphere_d<Traits> Min_sphere; typedef K::Point_d Point; const int n = 10; // number of points const int d = 5; // dimension of points int main () { Point P[n]; // n points double coord[d]; // d coordinates CGAL::Random r; // random number generator for (int i=0; i<n; ++i) { for (int j=0; j<d; ++j) coord[j] = r.get_double(); P[i] = Point(d, coord, coord+d); // random point } Min_sphere ms (P, P+n); // smallest enclosing sphere CGAL::set_pretty_mode (std::cout); std::cout << ms; // output the sphere return 0; }