Navigation: Up, Table of Contents, Bibliography, Index, Title Page

Definition

Let S(w) be a set of weighted points in 3. Let p(w)=(p,wp), p 3, wp and z(w)=(z,wz), z 3, wz be two weighted points. A weighted point p(w)=(p,wp) can also be seen as a sphere of center p and radius wp. The power product between p(w) and z(w) is defined as

(p(w),z(w)) = p-z 2-wp2-wz2

where p-z is the Euclidean distance between p and z. p(w) and z(w) are said to be orthogonal if (p(w)-z(w)) = 0 (see Figure reference).

Figure:  Orthogonal weighted points (picture in 2D)

Orthogonal weighted
points (picture in 2D)

Four weighted points have a unique common orthogonal weighted point called the power sphere. A sphere z(w) is said to be regular if p(w) S(w), (p(w)-z(w)) 0.

A triangulation of S(w) is regular if the power spheres of all simplices are regular.

The regular triangulation of S(w) is in fact the projection onto 3 of the convex hull of the four-dimensional points (p, p-O 2-wp), for p(w)=(p,wp) S(w). Note that all points of S(w) do not necessarily appear as vertices of the regular triangulation. To know more about regular triangulations, see for example [ES96].

When all weights are 0, power spheres are nothing more than circumscribing spheres, and the regular triangulation is exactly the Delaunay triangulation.

We will call the weighted point orthogonal to three weighted points in the plane defined by these three points the power circle. The power segment will denote the weighted point orthogonal to two weighted points on the line defined by these two points.

To simplify notation, p will often denote in the sequel either the point p 3 or the weighted point p(w)=(p,wp).

#include <CGAL/Regular_triangulation_3.h>

Inherits From

Triangulation_3<Traits,Tds>

Types

Inherits the types of Triangulation_3<Traits,Tds>.

typedef Traits::Bare_point
Bare_point; The type for points p of weighted points p(w)=(p,wp)
typedef Traits::Weighted_point
Weighted_point;

Creation

Regular_triangulation_3<Traits,Tds> rt;
Creates an empty regular triangulation.


Regular_triangulation_3<Traits,Tds> rt ( Traits traits);
Creates an empty regular triangulation with traits class traits.


Regular_triangulation_3<Traits,Tds> rt ( dt1);
Copy constructor.

Modifiers

Insertion

The following methods, which already exist in triangulations, are overloaded to ensure the property that all power spheres are regular.

Vertex_handle rt.insert ( Weighted_point p)
Inserts weighted point p in the triangulation. If this insertion creates a vertex, this vertex is returned. Otherwise, this method returns NULL.

Vertex_handle rt.insert ( Weighted_point p, Cell_handle start)
Same as the previous method, start is used as a starting place for the search done within the insertion.

The following method allows one to insert several points. It returns the number of inserted points.

template < class InputIterator >
int rt.insert ( InputIterator first, InputIterator last)
Inserts the points in the range [.first, last.).
Precondition: The value_type of first and last is Point.

Queries

Let us remark that

(p(w)-z(w)) > 0

is equivalent to

p lies outside the sphere with center z and radius sqrt(wp2+wz2).

This remark helps provide an intuition about the following predicates.

Figure:  side_of_power_circle

side\_of\_power\_circle

Bounded_side
rt.side_of_power_sphere ( Cell_handle c,
Weighted_point p)
Returns the position of the weighted point p with respect to the power sphere of c. More precisely, it returns:
- ON_BOUNDED_SIDE if (p(w)-z(c)(w))<0 where z(c)(w) is the power sphere of c. For an infinite cell this means either that p lies strictly in the half space limited by its finite facet and not containing any other point of the triangulation, or that the angle between p and the power circle of the finite facet of c is greater than /2.
- ON_BOUNDARY if p is orthogonal to the power sphere of c i.e. (p(w)-z(c)(w))=0. For an infinite cell this means that p is orthogonal to the power circle of its finite facet.
- ON_UNBOUNDED_SIDE if (p(w)-z(c)(w))>0 i.e. the angle between the weighted point p and the power sphere of c is less than /2 or if these two spheres do not intersect. For an infinite cell this means that p does not satisfy either of the two previous conditions.
Precondition: rt.dimension() =3.

Bounded_side rt.side_of_power_circle ( Facet f, Weighted_point p)
Returns the position of the point p with respect to the power circle of f. More precisely, it returns:
- in dimension 3:
- For a finite facet,
ON_BOUNDARY if p is orthogonal to the power circle in the plane of the facet,
ON_UNBOUNDED_SIDE when their angle is less than /2,
ON_BOUNDED_SIDE when it is greater than /2 (see Figure reference).
- For an infinite facet, it considers the plane defined by the finite facet of the cell f.first, and does the same as in dimension 2 in this plane.
- in dimension 2:
- For a finite facet,
ON_BOUNDARY if p is orthogonal to the circle,
ON_UNBOUNDED_SIDE when the angle between p and the power circle of f is less than /2, ON_BOUNDED_SIDE when it is greater than /2.
- For an infinite facet,
ON_BOUNDED_SIDE for a point in the open half plane defined by f and not containing any other point of the triangulation,
ON_UNBOUNDED_SIDE in the other open half plane.
If the point p is collinear with the finite edge e of f, it returns:
ON_BOUNDED_SIDE if (p(w)-z(e)(w))<0, where z(e)(w) is the power segment of e in the line supporting e,
ON_BOUNDARY if (p(w)-z(e)(w))=0,
ON_UNBOUNDED_SIDE if (p(w)-z(e)(w))>0 .
Precondition: rt.dimension() 2.

Bounded_side
rt.side_of_power_circle ( Cell_handle c,
int i,
Weighted_point p)
Same as the previous method for facet i of cell c.

Bounded_side
rt.side_of_power_segment ( Cell_handle c,
Weighted_point p)
In dimension 1, returns
ON_BOUNDED_SIDE if (p(w)-z(c)(w))<0, where z(c)(w) is the power segment of the edge represented by c,
ON_BOUNDARY if (p(w)-z(c)(w))=0,
ON_UNBOUNDED_SIDE if (p(w)-z(c)(w))>0 .
Precondition: rt.dimension() = 1.


begin of advanced section

Checking

bool rt.is_valid ( bool verbose = false)
Checks the combinatorial validity of the triangulation and the validity of its geometric embedding (see Section reference). Also checks that all the power spheres (resp. power circles in dimension 2, power segments in dimension 1) of cells (resp. facets in dimension 2, edges in dimension 1) are regular. When verbose is set to true, messages describing the first invalidity encountered are printed.
This method is mainly a debugging help for the users of advanced features.

end of advanced section


Next: Class declaration of Triangulation_Traits_3
Navigation: Up, Table of Contents, Bibliography, Index, Title Page
The GALIA project. Jan 18, 2000.