An object of the class Min_circle_2<Traits> is the unique circle of smallest area enclosing a finite (multi)set of points in two-dimensional Euclidean space 2. For a point set P we denote by mc(P) the smallest circle that contains all points of P. Note that mc(P) can be degenerate, i.e. mc(P)=∅ if P=∅ and mc(P)={p} if P={p}.
An inclusion-minimal subset S of P with mc(S)=mc(P) is called a support set, the points in S are the support points. A support set has size at most three, and all its points lie on the boundary of mc(P). In general, neither the support set nor its size are 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 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_circle_2<Traits> even if used only for points in two dimensions 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_circle_2<Traits> is not tuned for floating-point computations. The only advantage of Min_circle_2<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_circle_2<Traits> might still be an option in case your input number type cannot (efficiently) divide.
#include <CGAL/Min_circle_2.h>
The template parameter Traits is a model for MinCircle2Traits.
We provide the model CGAL::Min_circle_2_traits_2 using the two-dimensional Cgal kernel.
A Min_circle_2<Traits> object can be created from an arbitrary point set P and by specialized construction methods expecting no, one, two or three points as arguments. The latter methods can be useful for reconstructing mc(P) from a given support set S of P.
| |||
| |||
initializes min_circle to mc(P) with P being the set of points
in the range [first,last). If randomize is
true, a random permutation of P is computed in
advance, using the random numbers generator random.
Usually, this will not be necessary, however, the algorithm's
efficiency depends on the order in which the points are
processed, and a bad order might lead to extremely poor
performance (see example below).
| |||
| |||
initializes min_circle to
mc(∅), the empty set.
| |||
| |||
initializes min_circle to mc({p}), the set {p}.
| |||
| |||
initializes min_circle to mc({p1,p2}), the circle with diameter
equal to the segment connecting p1 and p2.
| |||
| |||
initializes min_circle to mc({p1,p2,p3}).
|
|
| |||
returns the number of points of min_circle, i.e. |P|. | ||||
|
| |||
returns the number of support points of min_circle, i.e. |S|. | ||||
|
| returns an iterator referring to the first point of min_circle. | ||
|
| returns the corresponding past-the-end iterator. | ||
|
| |||
returns an iterator referring to the first support point of min_circle. | ||||
|
| |||
returns the corresponding past-the-end iterator. | ||||
|
| |||
returns the i-th support point of min_circle. Between two
modifying operations (see below) any call to
min_circle.support_point(i) with the same i returns
the same point.
| ||||
|
| returns the current circle of min_circle. |
By definition, an empty Min_circle_2<Traits> has no boundary and no bounded side, i.e. its unbounded side equals the whole space 2.
|
| |
returns CGAL::ON_BOUNDED_SIDE, CGAL::ON_BOUNDARY, or CGAL::ON_UNBOUNDED_SIDE iff p lies properly inside, on the boundary of, or properly outside of min_circle, resp. | ||
|
| |
returns true, iff p lies properly inside min_circle. | ||
|
| |
returns true, iff p lies on the boundary of min_circle. | ||
|
| |
returns true, iff p lies properly outside of min_circle. | ||
|
| returns true, iff min_circle is empty (this implies degeneracy). |
|
| returns true, iff min_circle is degenerate, i.e. if min_circle is empty or equal to a single point, equivalently if the number of support points is less than 2. |
New points can be added to an existing min_circle, allowing to build mc(P) incrementally, e.g. if P is not known in advance. Compared to the direct creation of mc(P), this is not much slower, because the construction method is incremental itself.
|
| inserts p into min_circle and recomputes the smallest enclosing circle. | ||
| ||||
|
| |||
inserts the points in the range [first,last)
into min_circle and recomputes the smallest enclosing circle by
calling insert(p) for each point p in
[first,last).
| ||||
|
|
deletes all points in min_circle and sets min_circle to the empty set.
|
An object min_circle is valid, iff
|
| |
returns true, iff min_circle 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_circle to output stream os.
|
|
|
reads min_circle from input stream is.
|
CGAL::Min_ellipse_2<Traits>
CGAL::Min_sphere_d<Traits>
CGAL::Min_sphere_of_spheres_d<Traits>
CGAL::Min_circle_2_traits_2<K>
MinCircle2Traits
We implement the incremental algorithm of Welzl, with move-to-front heuristic [Wel91]. The whole implementation is described in [GS98a].
If randomization is chosen, 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 circle from scratch. The clear operation and the check for validity each takes linear time.
To illustrate the creation of Min_circle_2<Traits> and to show that randomization can be useful in certain cases, we give an example.
File: examples/Min_circle_2/min_circle_2.cpp
#include <CGAL/Exact_predicates_exact_constructions_kernel.h> #include <CGAL/Min_circle_2.h> #include <CGAL/Min_circle_2_traits_2.h> #include <iostream> // typedefs typedef CGAL::Exact_predicates_exact_constructions_kernel K; typedef CGAL::Min_circle_2_traits_2<K> Traits; typedef CGAL::Min_circle_2<Traits> Min_circle; typedef K::Point_2 Point; int main( int, char**) { const int n = 100; Point P[n]; for ( int i = 0; i < n; ++i) P[ i] = Point( (i%2 == 0 ? i : -i), 0); // (0,0), (-1,0), (2,0), (-3,0), ... Min_circle mc1( P, P+n, false); // very slow Min_circle mc2( P, P+n, true); // fast CGAL::set_pretty_mode( std::cout); std::cout << mc2; return 0; }