CGAL::Min_sphere_of_spheres_d<Traits>

Definition

An object of the class Min_sphere_of_spheres_d<Traits> is a data structure that represents the unique sphere of smallest volume enclosing a finite set of spheres in d-dimensional Euclidean space d. For a set S of spheres we denote by ms(S) the smallest sphere that contains all spheres of S; we call ms(S) the minsphere of S. ms(S) can be degenerate, i.e., ms(S)=∅, if S=∅ and ms(S)={s}, if S={s}. Any sphere in S may be degenerate, too, i.e., any sphere from S may be a point. Also, S may contain several copies of the same sphere.

An inclusion-minimal subset R of S with ms(R)=ms(S) is called a support set for ms(S); the spheres in R are the support spheres. A support set has size at most d+1, and all its spheres lie on the boundary of ms(S). (A sphere s' is said to lie on the boundary of a sphere s, if s' is contained in s and if their boundaries intersect.) In general, the support set is not unique.

The algorithm computes the center and the radius of ms(S), and finds a support set R (which remains fixed until the next insert(), clear() or set() operation). We also provide a specialization of the algorithm for the case when the center coordinates and radii of the input spheres are floating-point numbers. This specialized algorithm uses floating-point arithmetic only, is very fast and especially tuned for stability and robustness. Still, it's output may be incorrect in some (rare) cases; termination is guaranteed.

When default constructed, an instance of type Min_sphere_of_spheres_d<Traits> represents the set S=∅, together with its minsphere ms(S)=∅. You can add spheres to the set S by calling insert(). Querying the minsphere is done by calling the routines is_empty(), radius() and center_cartesian_begin(), among others.

In general, the radius and the Euclidean center coordinates of ms(S) need not be rational. Consequently, the algorithm computing the exact minsphere will have to deal with algebraic numbers. Fortunately, both the radius and the coordinates of the minsphere are numbers of the form ai+bit, where ai,bi,t and where t 0 is the same for all coordinates and the radius. Thus, the exact minsphere can be described by the number t, which is called the sphere's discriminant, and by d+1 pairs (ai,bi) 2 (one for the radius and d for the center coordinates).

#include <CGAL/Min_sphere_of_spheres_d.h>

Note: This class (almost) replaces CGAL::Min_sphere_d<Traits>, which solves the less general problem of finding the smallest enclosing ball of a set of points. Min_sphere_of_spheres_d<Traits> is faster than CGAL::Min_sphere_d<Traits>, and in contrast to the latter provides a specialized implementation for floating-point arithmetic which ensures correct results in a large number of cases (including highly degenerate ones). The only advantage of CGAL::Min_sphere_d<Traits> over 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, CGAL::Min_sphere_d<Traits> might still be an option in case your input number type cannot (efficiently) divide.

Requirements

The class Min_sphere_of_spheres_d<Traits> expects a model of the concept MinSphereOfSpheresTraits as its template argument.

Types

Min_sphere_of_spheres_d<Traits>::Sphere
is a typedef to Traits::Sphere.


Min_sphere_of_spheres_d<Traits>::FT
is a typedef to Traits::FT.


Min_sphere_of_spheres_d<Traits>::Result
is the type of the radius and of the center coordinates of the computed minsphere: When FT is an inexact number type (double, for instance), then Result is simply FT. However, when FT is an exact number type, then Result is a typedef to a derived class of std::pair<FT,FT>; an instance of this type represents the number a+b√t, where a is the first and b the second element of the pair and where the number t is accessed using the member function disciminant() of class Min_sphere_of_spheres_d<Traits>.


Min_sphere_of_spheres_d<Traits>::Algorithm
is either CGAL::LP_algorithm or CGAL::Farthest_first_heuristic. As is described in the documentation of concept MinSphereOfSpheresTraits, the type Algorithm reflects the method which is used to compute the minsphere. (Normally, Algorithm coincides with Traits::Algorithm. However, if the method Traits::Algorithm should not be supported anymore in a future release, then Algorithm will have another type.)


Min_sphere_of_spheres_d<Traits>::Support_iterator
non-mutable model of the STL concept BidirectionalIterator with value type Sphere. Used to access the support spheres defining the smallest enclosing sphere.


Min_sphere_of_spheres_d<Traits>::Cartesian_const_iterator
non-mutable model of the STL concept BidirectionalIterator to access the center coordinates of the minsphere.

Creation

Min_sphere_of_spheres_d<Traits> minsphere ( Traits traits = Traits());
creates a variable of type Min_sphere_of_spheres_d<Traits> and initializes it to ms(∅). If the traits parameter is not supplied, the class Traits must provide a default constructor.


template < typename InputIterator >
Min_sphere_of_spheres_d<Traits> minsphere ( InputIterator first, InputIterator last, Traits traits = Traits());
creates a variable minsphere of type Min_sphere_of_spheres_d<Traits> and inserts (cf. insert()) the spheres from the range [first,last).
Requirement: The value type of first and last is Sphere. If the traits parameter is not supplied, the class Traits must provide a default constructor.

Access Functions

Support_iterator minsphere.support_begin () returns an iterator referring to the first support sphere of minsphere.

Support_iterator minsphere.support_end () returns the corresponding past-the-end iterator.

FT minsphere.discriminant () returns the discriminant of minsphere. This number is undefined when FT is an inexact number type. When FT is exact, the center coordinates and the radius of the minsphere are numbers of the form a+b√t, where t is the discriminant of the minsphere as returned by this function.
Precondition: minsphere is not empty, and FT is an exact number type.

Result minsphere.radius () returns the radius of minsphere. If FT is an exact number type then the radius of the minsphere is the real number a+b√t, where t is the minsphere's discriminant, a is the first and b the second component of the pair returned by radius().
Precondition: minsphere is not empty.

Cartesian_const_iterator minsphere.center_cartesian_begin ()
returns a const-iterator to the first of the Traits::D center coordinates of minsphere. The iterator returns objects of type Result. If FT is an exact number type, then a center coordinate is represented by a pair (a,b) describing the real number a+b√t, where t is the minsphere's discriminant (cf. discriminant()).
Precondition: minsphere is not empty.

Cartesian_const_iterator minsphere.center_cartesian_end () returns the corresponding past-the-end iterator, i.e. center_cartesian_begin()+Traits::D.
Precondition: minsphere is not empty.

Predicates

bool minsphere.is_empty () returns true, iff minsphere is empty, i.e. iff ms(S)=∅.

Modifiers

void minsphere.clear () resets minsphere to ms(∅), with S:= ∅.

template < class InputIterator >
void minsphere.set ( InputIterator first, InputIterator last)
sets minsphere to the ms(S), where S is the set of spheres in the range [first,last).
Requirement: The value type of first and last is Sphere.

void minsphere.insert ( Sphere s) inserts the sphere s into the set S of instance minsphere.

template < class InputIterator >
void minsphere.insert ( InputIterator first, InputIterator last)
inserts the spheres in the range [first,last) into the set S of instance minsphere.
Requirement: The value type of first and last is Sphere.

Validity Check

An object minsphere is valid, iff

bool minsphere.is_valid () returns true, iff minsphere is valid. When FT is inexact, this routine always returns true.

Miscellaneous

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

See Also

CGAL::Min_sphere_d<Traits>
CGAL::Min_circle_2<Traits>

Implementation

We implement two algorithms, the LP-algorithm and a heuristic [MSW92]. As described in the documentation of concept MinSphereOfSpheresTraits, each has its advantages and disadvantages: Our implementation of the LP-algorithm has maximal expected running time O(2d n), while the heuristic comes without any complexity guarantee. In particular, the LP-algorithm runs in linear time for fixed dimension d. (These running times hold for the arithmetic model, so they count the number of operations on the number type FT.)

On the other hand, the LP-algorithm is, for inexact number types FT, much worse at handling degeneracies and should therefore not be used in such a case. (For exact number types FT, both methods handle all kinds of degeneracies.)

Currently, we require Traits::FT to be either an exact number type or double or float; other inexact number types are not supported at this time. Also, the current implementation only handles spheres with Cartesian coordinates; homogenous representation is not supported yet.

Example

File: examples/Min_sphere_of_spheres_d/min_sphere_of_spheres_d_d.cpp
#include <CGAL/Cartesian_d.h>
#include <CGAL/Random.h>
#include <CGAL/Gmpq.h>
#include <CGAL/Min_sphere_of_spheres_d.h>
#include <vector>

const int N = 1000;                       // number of spheres
const int D = 3;                          // dimension of points
const int LOW = 0, HIGH = 10000;          // range of coordinates and radii

typedef CGAL::Gmpq                        FT;
//typedef double                          FT;
typedef CGAL::Cartesian_d<FT>             K;
typedef CGAL::Min_sphere_of_spheres_d_traits_d<K,FT,D> Traits;
typedef CGAL::Min_sphere_of_spheres_d<Traits> Min_sphere;
typedef K::Point_d                        Point;
typedef Traits::Sphere                    Sphere;

int main () {
  std::vector<Sphere> S;                  // n spheres
  FT 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_int(LOW,HIGH);
    Point p(D,coord,coord+D);             // random center...
    S.push_back(Sphere(p,r.get_int(LOW,HIGH))); // ...and random radius
  }

  Min_sphere ms(S.begin(),S.end());       // check in the spheres
  CGAL_assertion(ms.is_valid());
}