Voronoi Tessellation output

In this section, we discuss the output given from voronoi. We consider a simple clipped example for examining this output.

using DelaunayTriangulation
using StableRNGs
rng = StableRNG(123)
p1 = (0.0, 1.0)
p2 = (3.0, -1.0)
p3 = (2.0, 0.0)
p4 = (-1.0, 2.0)
p5 = (4.0, 2.0)
p6 = (-2.0, -1.0)
p7 = (2.0, 1.0)
p8 = (5.0, 1.0)
points = [p1, p2, p3, p4, p5, p6, p7, p8]
tri = triangulate(points; rng)
vorn = voronoi(tri, clip = true, rng = rng)
vorn
Voronoi Tessellation.
    Number of generators: 8
    Number of polygon vertices: 21
    Number of polygons: 8

Now let's inspect vorn. The fields in vorn are:

propertynames(vorn)
(:triangulation, :generators, :polygon_points, :polygons, :circumcenter_to_triangle, :triangle_to_circumcenter, :unbounded_polygons, :cocircular_circumcenters, :adjacent, :boundary_polygons)

Let's examine this fields one at a time.

get_triangulation(vorn)

This field stores the triangulation object used to generate the Voronoi tessellation.

DelaunayTriangulation.get_triangulation(vorn)
Delaunay Triangulation.
   Number of vertices: 8
   Number of triangles: 9
   Number of edges: 16
   Has boundary nodes: false
   Has ghost triangles: true
   Curve-bounded: false
   Weighted: false
   Constrained: false

get_generators(vorn)

This field stores the generators of the Voronoi tessellation, i.e. the points associated with the polygons in the tessellation. These are simply the points in the triangulation, but are stored differently in case some points are not present in the triangulation.

DelaunayTriangulation.get_generators(vorn)
Dict{Int64, Tuple{Float64, Float64}} with 8 entries:
  5 => (4.0, 2.0)
  4 => (-1.0, 2.0)
  6 => (-2.0, -1.0)
  7 => (2.0, 1.0)
  2 => (3.0, -1.0)
  8 => (5.0, 1.0)
  3 => (2.0, 0.0)
  1 => (0.0, 1.0)

See that the generators are stored as a Dict, with the vertices mapping to their associated coordinates. The preferred way to access the generators is through get_generator, which can be used to obtain the coordinates for a given vertex. For example:

get_generator(vorn, 1)
(0.0, 1.0)
get_generator(vorn, 2, 6)
((3.0, -1.0), (-2.0, -1.0))

If you want to iterate over each generator, you should use each_generator, which simply returns an unordered iterator over the generators.

each_generator(vorn)
KeySet for a Dict{Int64, Tuple{Float64, Float64}} with 8 entries. Keys:
  5
  4
  6
  7
  2
  8
  3
  1

You can get the number of generators using num_generators:

DelaunayTriangulation.num_generators(vorn)
8

get_polygon_points(vorn)

This field is used for storing the coordinates of the Voronoi polygons:

DelaunayTriangulation.get_polygon_points(vorn)
21-element Vector{Tuple{Float64, Float64}}:
 (3.5, 0.5)
 (0.16666666666666666, -1.1666666666666665)
 (-1.5, 0.5)
 (1.0, 0.5)
 (0.5, -2.5)
 (1.0, 3.0)
 (1.5, 4.5)
 (4.0, 0.0)
 (4.5, 1.5)
 (5.0, 1.0)
 ⋮
 (2.0, -1.0)
 (3.0, -1.0)
 (-1.5, 0.5)
 (5.551115123125783e-17, 2.0)
 (-1.0, 2.0)
 (1.0, 2.0)
 (1.1102230246251565e-16, -1.0)
 (-2.0, -1.0)
 (0.24999999999999978, -1.0)

These coordinates are used to define the vertices of these polygons. It is important to note that these points are not guaranteed to be unique if a circumcenter appears on the boundary, or if the tessellation is clipped. To access the coordinate of a polygon associated with a given polygon vertex, use get_polygon_point:

get_polygon_point(vorn, 1)
(3.5, 0.5)

If you wanted the total number of polygon vertices, possibly counting duplicate vertices, you can use num_polygon_vertices:

num_polygon_vertices(vorn)
21

get_polygons(vorn)

This field is used for mapping polygon indices to their vertices, with these vertices referring to coordinates from get_polygon_points(vorn):

DelaunayTriangulation.get_polygons(vorn)
Dict{Int64, Vector{Int64}} with 8 entries:
  5 => [1, 9, 12, 11, 1]
  4 => [3, 16, 17, 3]
  6 => [20, 19, 3, 20]
  7 => [4, 1, 11, 18, 4]
  2 => [13, 14, 8, 1, 13]
  8 => [1, 8, 10, 9, 1]
  3 => [21, 13, 1, 4, 21]
  1 => [19, 21, 4, 18, 16, 3, 19]

For example, using get_polygon with the first polygon:

get_polygon(vorn, 1)
7-element Vector{Int64}:
 19
 21
  4
 18
 16
  3
 19

means that the first polygon has vertices [19, 21, 4, 18, 16, 3, 19]. All polygons are stored as a circular vector of counter-clockwise oriented vertices. Note also that polygon indices are the same as the generator vertex, so that e.g. this polygon 1 is that associated with generator 1. You can get the number of polygons using num_polygons:

num_polygons(vorn)
8

Of course, this is the same as the number of generators. If you wanted to get the coordinates instead of the vertices associated with a polygon, you use get_polygon_coordinates:

get_polygon_coordinates(vorn, 1)
7-element Vector{Tuple{Float64, Float64}}:
 (1.1102230246251565e-16, -1.0)
 (0.24999999999999978, -1.0)
 (1.0, 0.5)
 (1.0, 2.0)
 (5.551115123125783e-17, 2.0)
 (-1.5, 0.5)
 (1.1102230246251565e-16, -1.0)

This is almost the same result as doing

get_polygon_point.(Ref(vorn), get_polygon(vorn, 1))
7-element Vector{Tuple{Float64, Float64}}:
 (1.1102230246251565e-16, -1.0)
 (0.24999999999999978, -1.0)
 (1.0, 0.5)
 (1.0, 2.0)
 (5.551115123125783e-17, 2.0)
 (-1.5, 0.5)
 (1.1102230246251565e-16, -1.0)

but is more efficient.

When the tessellation is unbounded, this field is slightly different. Let's consider an example.

rng2 = StableRNG(321)
points2 = rand(rng2, 2, 10)
tri2 = triangulate(points2; rng = rng2)
vorn2 = voronoi(tri2, rng = rng2)
DelaunayTriangulation.get_polygons(vorn2)
Dict{Int64, Vector{Int64}} with 10 entries:
  5  => [12, 7, 10, -1, -3, 5, 12]
  4  => [8, 9, 3, 1, 4, 2, 8]
  6  => [5, 6, 8, 2, 12, 5]
  7  => [9, 13, -4, -5, 3, 9]
  2  => [8, 6, 13, 9, 8]
  10 => [10, 11, -2, -1, 10]
  9  => [-3, -4, 13, 6, 5, -3]
  8  => [7, 4, 1, 11, 10, 7]
  3  => [3, -5, -2, 11, 1, 3]
  1  => [12, 2, 4, 7, 12]

Notice that some of these polygons have negative vertices. These vertices are not the same as ghost vertices from triangulations. They do still represent points out at infinity, but their numbering is of no importance - only the fact that they are negative is. All polygons that have a negative vertex are unbounded, and the negative vertex is the vertex that is out at infinity. When we try to get the coordinates of a polygon with a negative vertex, we get an error:

julia> get_polygon_coordinates(vorn2, 5)ERROR: ArgumentError: The polygon is unbounded, so a bounding box must be provided. See DelaunayTriangulation.polygon_bounds for a reasonable default.

We need to provide a bounding box to obtain these coordinates. As mentioned, we can obtain a reasonable default for a bounding box using polygon_bounds, which returns the bounding box in the form (xmin, xmax, ymin, ymax):

bbox = polygon_bounds(vorn2)
(-0.0845250117400761, 0.7292550055392539, -0.07703000058813853, 0.8087847404970718)

If we now use this bounding box to get the coordinates of the polygon, we get a polygon that will have been clipped with this box:

get_polygon_coordinates(vorn2, 4, bbox)
8-element Vector{Tuple{Float64, Float64}}:
 (0.3800182607228671, 0.4195812181251387)
 (0.5248526830200166, 0.30356339617239647)
 (0.7292550055392539, 0.41037087429841546)
 (0.7292550055392539, 0.41037087429841546)
 (0.6204044502699297, 0.543145666748361)
 (0.5062082223280505, 0.6112493082368492)
 (0.36237000338690506, 0.5169933977507611)
 (0.3800182607228671, 0.4195812181251387)

For iterating over polygons, there are two main methods. If you only care about the vertices and not the index associated with a polygon, you can use each_polygon:

each_polygon(vorn)
ValueIterator for a Dict{Int64, Vector{Int64}} with 8 entries. Values:
  [1, 9, 12, 11, 1]
  [3, 16, 17, 3]
  [20, 19, 3, 20]
  [4, 1, 11, 18, 4]
  [13, 14, 8, 1, 13]
  [1, 8, 10, 9, 1]
  [21, 13, 1, 4, 21]
  [19, 21, 4, 18, 16, 3, 19]

Alternatively, each_polygon_index can be used to iterate over the index of each polygon (equivalently, over each generator vertex), followed by get_polygon to access the associated vertices of that polygon:

for v in each_polygon_index(vorn)
    println(v, ": ", get_polygon(vorn, v))
end
5: [1, 9, 12, 11, 1]
4: [3, 16, 17, 3]
6: [20, 19, 3, 20]
7: [4, 1, 11, 18, 4]
2: [13, 14, 8, 1, 13]
8: [1, 8, 10, 9, 1]
3: [21, 13, 1, 4, 21]
1: [19, 21, 4, 18, 16, 3, 19]

get_circumcenter_to_triangle(vorn)

When a Voronoi tessellation is constructed, the circumcenters of the triangles from the underlying Delaunay triangulation are used to define the polygon vertices. This field is used to associate a given polygon vertex with the triangle that it came from:

DelaunayTriangulation.get_circumcenter_to_triangle(vorn)
Dict{Int64, Tuple{Int64, Int64, Int64}} with 12 entries:
  5  => (3, 6, 2)
  7  => (7, 5, 4)
  -3 => (8, 2, -1)
  1  => (8, 5, 3)
  4  => (3, 7, 1)
  6  => (7, 4, 1)
  -5 => (5, 8, -1)
  -1 => (2, 6, -1)
  2  => (6, 3, 1)
  -2 => (6, 4, -1)
  -4 => (4, 5, -1)
  3  => (4, 6, 1)

We see, for example, that the fifth polygon vertex is derived from the circumcenter of the triangle (3, 6, 2). We can verify this:

p = DelaunayTriangulation.triangle_circumcenter(tri, (3, 6, 2))
q = get_polygon_point(vorn, 5)
p == q
true

The negative vertices show what ghost triangle is being associated with an unbounded edge of the tessellation. Since we clipped our tessellation, these no longer appear in the previously mentioned fields, but they are still stored in this field; they are used in clipping for determining the direction of the unbounded edges. If you wanted the triangle associated with a specific vertex, you could use

DelaunayTriangulation.get_circumcenter_to_triangle(vorn, 5)
(3, 6, 2)

get_triangle_to_circumcenter(vorn)

This map is the inverse of get_circumcenter_to_triangle(vorn), and is used to map a triangle to its circumcenter:

DelaunayTriangulation.get_triangle_to_circumcenter(vorn)
Dict{Tuple{Int64, Int64, Int64}, Int64} with 14 entries:
  (8, 5, 3)  => 1
  (2, 6, -1) => -1
  (6, 4, -1) => -2
  (3, 7, 1)  => 4
  (5, 8, -1) => -5
  (6, 3, 1)  => 2
  (4, 6, 1)  => 3
  (3, 6, 2)  => 5
  (7, 5, 4)  => 7
  (5, 7, 3)  => 1
  (8, 2, -1) => -3
  (4, 5, -1) => -4
  (8, 3, 2)  => 1
  (7, 4, 1)  => 6

The triangles are stored so that the minimum vertex is always last, maintaining the counter-clockwise orientation of each. You can get the circumcenter index associated with a given triangle using:

DelaunayTriangulation.get_triangle_to_circumcenter(vorn, (3, 7, 1))
4

get_unbounded_polygons(vorn)

This field is used to store the polygons that are unbounded. For our clipped tessellation, this field is empty:

DelaunayTriangulation.get_unbounded_polygons(vorn)
Set{Int64}()

For our unclipped tessellation, we have

DelaunayTriangulation.get_unbounded_polygons(vorn2)
Set{Int64} with 5 elements:
  5
  7
  10
  9
  3

and we can verify that this is consistent with our observation that each unbounded polygon has a negative vertex:

S1 = DelaunayTriangulation.get_unbounded_polygons(vorn2)
S2 = Set{Int}()
for v in each_polygon_index(vorn2)
    C = get_polygon(vorn2, v)
    if any(<(0), C)
        push!(S2, v)
    end
end
S1 == S2
true

get_cocircular_circumcenters(vorn)

In cases where pairs of triangles are cocircular, meaning lie on the same circle, their circumcenter will be equal and thus the polygon vertex associated with these two triangles will be the same. This field will store the vertices of these circumcenters that came from cocircular triangles. For our case, there is a single entry:

DelaunayTriangulation.get_cocircular_circumcenters(vorn)
Set{Int64} with 1 element:
  1

This means that the first polygon vertex is associated with at least two cocircular triangles, and is why, for example, the values from get_triangle_to_circumcenter(vorn) are not unique:

t2c = DelaunayTriangulation.get_triangle_to_circumcenter(vorn)
allunique(values(t2c))
false

and also explains why get_circumcenter_to_triangle(vorn) and get_triangle_to_circumcenter(vorn) are not perfect inverses of each other. The triangle mapped to from 5 inside the circumcenter_to_triangle map is simply the first one that the algorithm happened to encounter. To find the two cocircular triangles, we could use the triangle_to_circumcenter map:

triangles = NTuple{3,Int}[]
for (t, c) in t2c
    c == 1 && push!(triangles, t)
end

In fact, there are three triangles that all have the same circumcenter!

DelaunayTriangulation.triangle_circumcenter.(Ref(tri), triangles)
3-element Vector{Tuple{Float64, Float64}}:
 (3.5, 0.5)
 (3.5, 0.5)
 (3.5, 0.5)

get_adjacent(vorn)

This map is similar to the Adjacent map for triangulations:

get_adjacent(vorn)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 33 entries:
  (18, 16) => 1
  (11, 18) => 7
  (19, 21) => 1
  (19, 3)  => 6
  (14, 8)  => 2
  (1, 9)   => 5
  (4, 21)  => 3
  (16, 17) => 4
  (13, 14) => 2
  (21, 13) => 3
  (21, 4)  => 1
  (9, 1)   => 8
  (10, 9)  => 8
  (3, 19)  => 1
  (17, 3)  => 4
  (1, 4)   => 3
  (1, 13)  => 2
  (9, 12)  => 5
  (4, 1)   => 7
  ⋮        => ⋮

In this case, oriented edges are being mapped to the polygon it belongs to. For example,

get_adjacent(vorn, 3, 19)
1

means that the edge (3, 19) belongs to the first polygon; these vertices 3 and 19 refer to the vertices from the polygon, while the 1 from the output refers to a generator.

get_boundary_polygons(vorn)

This field is used to store the polygons that are on the boundary of the tessellation, and is only relevant for clipped tessellations. For our clipped tessellation, we have:

DelaunayTriangulation.get_boundary_polygons(vorn)
Set{Int64} with 8 elements:
  5
  4
  6
  7
  2
  8
  3
  1

For our unclipped tessellation, the field is empty:

DelaunayTriangulation.get_boundary_polygons(vorn2)
Set{Int64}()