Triangulation Output

In this section, we discuss the output given from triangulate. Let's take a curve-bounded domain and inspect its output in detail; for information about this domain in particular, see the curve-bounded domain tutorial. First, here is the triangulation.

using DelaunayTriangulation
using StableRNGs
using DelaunayTriangulation: EllipticalArc
curve = [
    [
        [1, 2, 3], [EllipticalArc((2.0, 0.0), (-2.0, 0.0), (0.0, 0.0), 2, 1 / 2, 0.0)]
    ],
    [
        [BSpline([(0.0, 0.4), (1.0, 0.2), (0.0, 0.1), (-1.0, 0.2), (0.0, 0.4)])]
    ],
    [
        [4, 5, 6, 7, 4]
    ],
    [
        [BezierCurve([(0.0, -2.0), (0.0, -2.5), (-1.0, -2.5), (-1.0, -3.0)])], [CatmullRomSpline([(-1.0, -3.0), (0.0, -4.0), (1.0, -3.0), (0.0, -2.0)])]
    ],
    [
        [12, 11, 10, 12]
    ],
    [
        [CircularArc((1.1, -3.0), (1.1, -3.0), (0.0, -3.0), positive=false)]
    ]
]
points = [(-2.0, 0.0), (0.0, 0.0), (2.0, 0.0), (-2.0, -5.0), (2.0, -5.0), (2.0, -1 / 10), (-2.0, -1 / 10), (-1.0, -3.0), (0.0, -4.0), (0.0, -2.3), (-0.5, -3.5), (0.9, -3.0)]
rng = StableRNG(123)
tri = triangulate(points; boundary_nodes=curve, rng)
refine!(tri; max_area=1e-2, rng)
tri
Delaunay Triangulation.
   Number of vertices: 1920
   Number of triangles: 3433
   Number of edges: 5353
   Has boundary nodes: true
   Has ghost triangles: true
   Curve-bounded: true
   Weighted: false
   Constrained: true

Now let's inspect tri. If we look at the fields in tri, we see that there is a lot of information stored in tri:

propertynames(tri)
(:points, :triangles, :boundary_nodes, :interior_segments, :all_segments, :weights, :adjacent, :adjacent2vertex, :graph, :boundary_curves, :boundary_edge_map, :ghost_vertex_map, :ghost_vertex_ranges, :convex_hull, :representative_point_list, :polygon_hierarchy, :boundary_enricher, :cache)

Note that each field X has an associated accessor DelaunayTriangulation.get_X. Let's go through each field.

Geometry Fields

First, we list the fields relating to the actual geometry.

get_points(tri)

This field stores the points that were triangulated. This may be a mutated version (not a copy, though, tri.points === points above) of the points provided into tri, note.

get_points(tri)
2066-element Vector{Tuple{Float64, Float64}}:
 (-2.0, 0.0)
 (0.0, 0.0)
 (2.0, 0.0)
 (-2.0, -5.0)
 (2.0, -5.0)
 (2.0, -0.1)
 (-2.0, -0.1)
 (-1.0, -3.0)
 (0.0, -4.0)
 (0.0, -2.3)
 ⋮
 (-1.7471596889265435, -0.9910059470419887)
 (1.3662028137779731, -3.8247220533036566)
 (0.6691499520988589, -0.36911179341389805)
 (-0.8995276363058041, -1.4188380367630027)
 (0.4820683876129333, -4.7068525042954175)
 (1.698256524308673, -1.1015931034194562)
 (-0.8402653313886871, -0.4170053221048895)
 (1.539769650409285, -2.0379365189133525)
 (1.8033813809067119, -0.5255248910528076)

In some cases, the points in this vector will not all appear in tri, which is why it is recommend you work with the vertices instead, via each_vertex(tri) or, for the solid or ghost vertices use each_solid_vertex(tri) or each_ghost_vertex(tri), respectively.

get_triangles(tri)

This field stores all the triangles in the triangulation, including both solid and ghost triangles.

get_triangles(tri)
Set{Tuple{Int64, Int64, Int64}} with 3840 elements:
  (445, 250, 215)
  (73, 36, -6)
  (1386, 1385, 1384)
  (1624, 736, 1622)
  (1817, 743, 850)
  (1717, 977, 979)
  (1971, 1064, 696)
  (1880, 1877, 798)
  (1477, 1474, 1468)
  (1945, 718, 1178)
  (1788, 729, 845)
  (2051, 1026, 1979)
  (1731, 1084, 1071)
  (2010, 787, 1706)
  (1379, 744, 1378)
  (294, 121, 268)
  (12, 558, -7)
  (443, 120, 218)
  (1341, 1338, 948)
  ⋮ 

This is also the same as each_triangle(tri). If you only wanted the solid triangles, you could use each_solid_triangle(tri), and similarly for the ghost triangles with each_ghost_triangle(tri). All the triangles in these sets are defined to be counter-clockwise:

using DelaunayTriangulation: is_positively_oriented, triangle_orientation
all(T -> is_positively_oriented(triangle_orientation(get_point(tri, triangle_vertices(T)...)...)), each_solid_triangle(tri))
true

The above check is not true for all the ghost triangles since the domain is non-convex.

get_boundary_nodes(tri)

This field stores the boundary of the triangulation.

get_boundary_nodes(tri)
6-element Vector{Vector{Vector{Int64}}}:
 [[1, 524, 59, 1661, 58, 606, 299, 1773, 57, 345  …  352, 54, 576, 298, 1570, 55, 1659, 56, 528, 3], [3, 17, 529, 48, 2026, 510, 1572, 60, 1578, 207  …  208, 1759, 50, 535, 504, 2028, 49, 525, 18, 1]]
 [[13, 1819, 308, 154, 98, 97, 423, 96, 95, 20  …  21, 86, 87, 432, 88, 89, 102, 109, 1818, 13]]
 [[4, 1259, 1255, 1260, 955, 1805, 1245, 1402, 523, 1758  …  1201, 69, 1226, 1025, 1250, 958, 1252, 1234, 1256, 4]]
 [[14, 23, 150, 478, 536, 540, 22, 1647, 515, 473, 471, 24, 8], [8, 40, 41, 42, 43, 44, 79, 80, 78, 77  …  556, 562, 91, 559, 549, 588, 38, 592, 369, 14]]
 [[12, 1952, 415, 174, 1708, 364, 65, 616, 305, 544  …  479, 10, 551, 196, 586, 66, 611, 67, 558, 12]]
 [[15, 907, 1748, 737, 1520, 968, 1525, 31, 1375, 922  …  1716, 34, 1719, 847, 1959, 731, 1593, 778, 1590, 15]]

Notice that these boundary nodes are not the same as those provided into triangulate above. Instead of having curve, triangulate constructs a new boundary_nodes vector and populates that with the piecewise linear approximation. You could work with these nodes using get_boundary_nodes again. For examples, the nodes associated with the first curve are given by

get_boundary_nodes(tri, 1)
2-element Vector{Vector{Int64}}:
 [1, 524, 59, 1661, 58, 606, 299, 1773, 57, 345  …  352, 54, 576, 298, 1570, 55, 1659, 56, 528, 3]
 [3, 17, 529, 48, 2026, 510, 1572, 60, 1578, 207  …  208, 1759, 50, 535, 504, 2028, 49, 525, 18, 1]

The second section of this curve is given by

get_boundary_nodes(tri, 1, 2)
39-element Vector{Int64}:
    3
   17
  529
   48
 2026
  510
 1572
   60
 1578
  207
    ⋮
 1759
   50
  535
  504
 2028
   49
  525
   18
    1

get_interior_segments(tri)

This field stores the interior segments of the triangulation. These are the segments that are not part of the boundary.

get_interior_segments(tri)
Set{Tuple{Int64, Int64}}()

In our case, there are no such segments, but there would be if we had used segments when calling triangulate. We could add a segment now, for example

r = DelaunayTriangulation.num_points(tri)
add_point!(tri, -3/2, -4.0, concavity_protection = true)
add_point!(tri, -3/2, -1.0, concavity_protection = true)
add_segment!(tri, r + 1, r + 2)
Delaunay Triangulation.
   Number of vertices: 1922
   Number of triangles: 3437
   Number of edges: 5359
   Has boundary nodes: true
   Has ghost triangles: true
   Curve-bounded: true
   Weighted: false
   Constrained: true

Now we see that we have some segments:

get_interior_segments(tri)
Set{Tuple{Int64, Int64}} with 1 element:
  (2067, 2068)

Note that, for any segment (i, j), only one of (i, j) or (j, i) will appear in the list of interior segments.

get_all_segments(tri)

In contrast to the interior_segments field, which only stores the interior segments, this field stores both the interior and boundary segments.

get_all_segments(tri)
Set{Tuple{Int64, Int64}} with 408 elements:
  (159, 479)
  (1634, 721)
  (1578, 207)
  (1245, 1402)
  (1854, 972)
  (155, 309)
  (117, 118)
  (92, 93)
  (1819, 308)
  (1315, 596)
  (4, 1259)
  (1035, 1723)
  (528, 3)
  (639, 1489)
  (36, 73)
  (1852, 1027)
  (103, 2)
  (562, 91)
  (1209, 714)
  ⋮ 

Just as for the interior segments, only one of (i, j) or (j, i) will appear in the list of all segments.

get_weights(tri)

This field stores the weights associated with each point in the triangulation. Currently, this field is not used properly anywhere and should not be inspected until weighted triangulations are fully implemented. By default, you will see a ZeroWeight():

get_weights(tri)
DelaunayTriangulation.ZeroWeight()

Topology Fields

Now we list the fields relating to the topology of the triangulation itself.

get_adjacent(tri)

This field stores the adjacent map, mapping each edge (u, v) in the triangulation to the vertex w such that (u, v, w) is a positively oriented triangle in the triangulation. In cases where there is no such triangle, the vertex 0 is returned. For this triangulation, we have:

get_adjacent(tri)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 11532 entries:
  (943, 1770)  => -8
  (1701, 984)  => 1145
  (71, 1598)   => 1691
  (1264, 1232) => 1235
  (975, 1734)  => 1673
  (306, 456)   => 459
  (953, 1333)  => 1329
  (1186, 1171) => 715
  (1079, 1416) => 1926
  (423, 97)    => -3
  (1949, 141)  => 416
  (231, 181)   => 444
  (1160, 1096) => 1070
  (1308, 1309) => 758
  (1597, 1147) => 1640
  (1454, 2021) => 1249
  (992, 1934)  => 1748
  (1558, 1473) => 1898
  (1317, 1041) => 1318
  ⋮            => ⋮

For example, the mapping (1648, 2064) => 803 implies that the triangle (1648, 2064, 803) is in the triangulation and is positively oriented, which we can verify:

DelaunayTriangulation.contains_triangle(tri, 1648, 2064, 803)
((1648, 2064, 803), false)

It is important to note that, for any triangle (u, v, w), the mappings (u, v) => w, (v, w) => u, and (w, u) => v will all appear in the adjacent map. For example:

get_adjacent(tri, 1648, 2064), get_adjacent(tri, 2064, 803), get_adjacent(tri, 803, 1648)
(0, 0, 0)

You can use get_adjacent(tri, u, v) to find this vertex w, as we demonstrated above. For cases where the edge does not exist, you will get 0:

get_adjacent(tri, 1, 2)
0

One last thing to notice is that some of the vertices w will be ghost vertices, and similarly some of the edges (u, v) will be ghost edges. For example,

get_adjacent(tri, 423, 97)
-3

This vertex -3 means that (423, 97) is an edge of the boundary associated with the ghost vertex -3.

get_adjacent2vertex(tri)

The adjacent2vertex field is similar to adjacent. The difference is that, instead of mapping edges to vertices, adjacent2vertex maps vertices w to the set of all edges (u, v) such that (u, v, w) is a positively oriented triangle in the triangulation. For this triangulation, we have:

get_adjacent2vertex(tri)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 1930 entries:
  1144 => Set([(754, 1704), (1143, 754), (755, 1145), (1145, 695), (695, 1143),…
  1175 => Set([(2002, 1059), (1059, 1620), (718, 1174), (1174, 2002), (1620, 18…
  1953 => Set([(956, 1281), (1936, 1274), (1872, 956), (1275, 1936), (1274, 187…
  719  => Set([(1632, 1656), (1055, 1687), (1687, 2056), (1162, 1055), (2056, 1…
  1546 => Set([(1744, 848), (883, 749), (848, 1811), (1545, 1744), (749, 1545),…
  1703 => Set([(685, 1306), (1704, 685), (755, 1704), (1702, 755), (1306, 756),…
  1956 => Set([(1646, 641), (1143, 1646), (754, 1143), (641, 1112), (1112, 754)…
  1028 => Set([(1185, 1187), (1891, 1183), (1182, 1185), (1189, 1891), (1187, 1…
  699  => Set([(981, 983), (852, 981), (983, 985), (986, 1779), (985, 986), (17…
  831  => Set([(1766, 832), (1806, 836), (832, 1961), (622, 1766), (1961, 1806)…
  1299 => Set([(1670, 1297), (1286, 1241), (1240, 1286), (1241, 1670), (1298, 1…
  1438 => Set([(1832, 1743), (1433, 1832), (1835, 1437), (1432, 1929), (1929, 1…
  -5   => Set([(8, 24), (150, 23), (22, 540), (478, 150), (536, 478), (1647, 22…
  1074 => Set([(1068, 1730), (1075, 1068), (1838, 1075), (2041, 1838), (1730, 2…
  319  => Set([(133, 127), (318, 133), (321, 320), (156, 321), (355, 156), (127…
  687  => Set([(1107, 1108), (1105, 1106), (1721, 1104), (1104, 1105), (1106, 1…
  1812 => Set([(688, 1813), (1044, 2044), (1042, 2061), (1813, 1042), (2044, 68…
  1199 => Set([(1021, 1203), (1203, 2036), (1200, 1021), (2036, 1208), (2048, 1…
  185  => Set([(271, 314), (442, 441), (268, 442), (441, 271), (317, 268), (314…
  ⋮    => ⋮

An example of this mapping is:

get_adjacent2vertex(tri, 719)
Set{Tuple{Int64, Int64}} with 6 elements:
  (1632, 1656)
  (1055, 1687)
  (1687, 2056)
  (1162, 1055)
  (2056, 1632)
  (1656, 1162)

This output means that (2057, 1625, 719), (1055, 1680, 719), (1625, 1649, 719), (1649, 1162, 719), (1162, 1055, 719), and (1680, 2057, 719) are all positively oriented triangles in the triangulation, and these are the only triangles that contain the vertex 719. In contrast to get_adjacent, calling get_adjacent2vertex on a vertex not in the triangulation will throw a KeyError. For ghost vertices, you will get the set of all edges on the boundary associated with that vertex, for example

get_adjacent2vertex(tri, -1)
Set{Tuple{Int64, Int64}} with 32 elements:
  (2, 103)
  (54, 352)
  (533, 99)
  (3, 528)
  (56, 1659)
  (298, 576)
  (299, 606)
  (568, 2)
  (138, 285)
  (61, 533)
  (99, 568)
  (1773, 299)
  (103, 240)
  (59, 524)
  (58, 1661)
  (290, 136)
  (524, 1)
  (576, 54)
  (1570, 298)
  ⋮ 

gives the set of all edges on the boundary associated with the ghost vertex -1. It is important to note that the edges in this set are not returned in any particular order, and so you should not rely on the order of the edges in the set.

get_graph(tri)

The last field relating to topology is the graph, which stores the graph information of the underlying triangulation. In particular, it stores the information about which vertices share an edge. All together, we have three different types of connectivity information being stored in a triangulation:

  1. Adjacent: The map from edges to vertices.
  2. Adjacent2Vertex: The map from vertices to edges.
  3. Graph: The map from vertices to vertices.

For this triangulation, we have:

get_graph(tri)
Graph
    Number of edges: 5770
    Number of vertices: 1930

This output by itself of course isn't too useful. The returned Graph does have its own fields,

propertynames(get_graph(tri))
(:vertices, :edges, :neighbours)

which you can expect using DelaunayTriangulation.get_vertices(tri), DelaunayTriangulation.get_edges(tri), and get_neighbours(tri). Let us go through each of these fields one at a time.

each_solid_vertex(tri)
EachSolidVertex iterator over vertices in a triangulation.
collect(each_solid_vertex(tri))
1922-element Vector{Int64}:
 1144
 1175
 1953
  719
 1546
 1703
 1956
 1028
  699
  831
    ⋮
  641
  979
 1711
   65
 1827
 2054
  298
 1402
 1115

Since these iterators are derived from Sets, you should not rely on the particular ordering of vertices returned.

each_edge(tri)
Set{Tuple{Int64, Int64}} with 5770 elements:
  (1964, 724)
  (2047, 2046)
  (1701, 984)
  (36, 73)
  (353, 348)
  (1527, 767)
  (71, 1598)
  (448, 187)
  (1264, 1232)
  (424, 333)
  (1186, 1171)
  (1949, 141)
  (231, 181)
  (1160, 1096)
  (1597, 1147)
  (575, 212)
  (461, 450)
  (187, 233)
  (1317, 1041)
  ⋮ 
collect(each_edge(tri))
5770-element Vector{Tuple{Int64, Int64}}:
 (1964, 724)
 (2047, 2046)
 (1701, 984)
 (36, 73)
 (353, 348)
 (1527, 767)
 (71, 1598)
 (448, 187)
 (1264, 1232)
 (424, 333)
 ⋮
 (1415, 678)
 (1770, 1769)
 (752, -4)
 (67, 558)
 (310, 13)
 (2044, 1812)
 (1398, 1243)
 (1715, 584)
 (192, 179)

Since these iterators are derived from Sets, you should not rely on the particular ordering of edges returned.

  • The neighbours are typically the more useful part of Graph. Looking at get_neighbours by itself, you obtain a mapping:
get_neighbours(tri)
Dict{Int64, Set{Int64}} with 1930 entries:
  1144 => Set([1704, 754, 695, 755, 1145, 1143])
  1175 => Set([1620, 1803, 1059, 1890, 1174, 2002, 718])
  1953 => Set([1275, 1936, 1872, 1274, 956, 1281])
  719  => Set([1656, 1162, 1632, 2056, 1055, 1687])
  1546 => Set([1744, 883, 848, 1811, 749, 1545])
  1703 => Set([1704, 755, 1306, 1702, 685, 756])
  1956 => Set([754, 641, 1112, 1143, 1646])
  1028 => Set([1189, 1891, 1182, 1183, 1185, 1187])
  699  => Set([985, 852, 983, 986, 981, 1779])
  831  => Set([836, 622, 832, 1961, 1806, 1766])
  1299 => Set([1241, 1240, 1297, 1298, 1670, 1286])
  1438 => Set([1437, 1433, 1832, 1743, 1835, 1929, 1432])
  -5   => Set([515, 24, 8, 1647, 23, 22, 473, 478, 14, 471, 536, 540, 150])
  1074 => Set([1068, 1838, 1075, 1730, 2041])
  319  => Set([318, 133, 321, 355, 127, 156, 320])
  687  => Set([1107, 1106, 1104, 1108, 1721, 1105])
  1812 => Set([1813, 2044, 2061, 1044, 688, 1042, 1043])
  1199 => Set([2036, 1203, 1208, 2048, 1200, 1021])
  185  => Set([441, 442, 268, 314, 271, 317, 315])
  ⋮    => ⋮

This mapping shows the neighbours for each vertex in the triangulation. For example,

get_neighbours(tri, 1)
Set{Int64} with 4 elements:
  -1
  18
  524
  -2

shows that the vertices of 1 are -1, 18, 524, and -2. These two ghost vertices, -1 and -2, imply that 1 is on the boundary of the triangulation and is also at the corner of the two sections of the boundary associated with the ghost vertices -1 and -2. Similarly,

get_neighbours(tri, -2)
Set{Int64} with 39 elements:
  16
  2028
  105
  60
  430
  410
  17
  1
  83
  49
  407
  510
  1578
  529
  338
  2026
  3
  51
  84
  ⋮ 

shows the set of all vertices that are on the boundary associated with the ghost vertex -2. Once again, the vertices in this set are not returned in any particular order, and so you should not rely on the order of the vertices in the set. If you did want to traverse the boundary in some order, you would need to use get_adjacent to do so, or use the boundary_nodes; see, for example, the get_triangulation_area function in this tutorial. (If you have a node on the boundary you are interested in, DelaunayTriangulation.get_left_boundary_node and DelaunayTriangulation.get_right_boundary_node are available.)

Boundary Handling Fields

The next set of fields relate to how the boundary is handled in the triangulation.

get_boundary_curves(tri)

This field is relevant only for curve-bounded domains, like the triangulation in this example. It stores the curves that were used to define the boundary of the triangulation. For our triangulation, we have

get_boundary_curves(tri)
(PiecewiseLinear(), EllipticalArc, BSpline, PiecewiseLinear(), BezierCurve, CatmullRomSpline, PiecewiseLinear(), CircularArc)

This output is not just giving the names of the curves - the actual curves themselves are stored. For example,

get_boundary_curves(tri)[3]
(::BSpline) (generic function with 1 method)

is the actual BSpline that we have provided. The curves are stored in the same order as they were provided, so that the jth curve is associated with the jth section (i.e., the ghost vertex -j) of the boundary. The parts of the boundary that are a sequence of straight lines between vertices are represented by a PiecewiseLinear curve. For a triangulation that is not curve-bounded, the output is an empty Tuple. For example:

tri2 = triangulate(rand(2, 50))
get_boundary_curves(tri2)
()

get_boundary_edge_map(tri)

This field is used to give information about where boundary edges are located in the triangulation. The precise definition is a bit cumbersome to describe: For a given boundary edge (u, v), the boundary_edge_map will map (u, v) to a tuple of the form (pos, ℓ), so that pos is the position of the edge in the boundary, and is the length of the boundary up to the edge. In particular, if get_boundary_edge_map(tri, u, v) == (pos, ℓ) and bn = get_boundary_nodes(tri, pos), then get_boundary_nodes(bn, ℓ) = u and get_boundary_nodes(bn, ℓ + 1) = v. For our triangulation, we have

get_boundary_edge_map(tri)
Dict{Tuple{Int64, Int64}, Tuple{Tuple{Int64, Int64}, Int64}} with 407 entries:
  (159, 479)   => ((5, 1), 23)
  (1634, 721)  => ((6, 1), 34)
  (1578, 207)  => ((1, 2), 9)
  (1245, 1402) => ((3, 1), 7)
  (1854, 972)  => ((3, 1), 30)
  (155, 309)   => ((1, 2), 18)
  (117, 118)   => ((4, 2), 33)
  (92, 93)     => ((4, 2), 27)
  (1819, 308)  => ((2, 1), 2)
  (1315, 596)  => ((3, 1), 86)
  (4, 1259)    => ((3, 1), 1)
  (1035, 1723) => ((6, 1), 41)
  (528, 3)     => ((1, 1), 32)
  (639, 1489)  => ((3, 1), 75)
  (36, 73)     => ((4, 2), 20)
  (1852, 1027) => ((6, 1), 32)
  (103, 2)     => ((1, 1), 16)
  (562, 91)    => ((4, 2), 69)
  (1209, 714)  => ((3, 1), 113)
  ⋮            => ⋮

Since our domain has multiple curves, the pos values are Tuples of the form (m, n), where m is the curve and n is the section of that curve. So, for example, we have

pos, ℓ = get_boundary_edge_map(tri, 1493, 763)

This means that (1493, 763) is the 45th edge of the first section on the third boundary:

bn = get_boundary_nodes(tri, pos) # notice that we can pass Tuples (m, n) as a single argument
get_boundary_nodes(bn, ℓ), get_boundary_nodes(bn, ℓ + 1)

For simpler domains, pos may be either a vector of vertices (in the case of a contiguous boundary) or a single vertex (in the case of a single boundary). For the contiguous case, we have:

tri3 = triangulate_rectangle(0, 1, 0, 1, 10, 10, single_boundary = true)
get_boundary_edge_map(tri3)
Dict{Tuple{Int64, Int64}, Tuple{Vector{Int64}, Int64}} with 36 entries:
  (91, 81)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (4, 5)    => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (1, 2)    => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (96, 95)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (60, 70)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (9, 10)   => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (98, 97)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (21, 11)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (100, 99) => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (95, 94)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (93, 92)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (30, 40)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (31, 21)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (61, 51)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (7, 8)    => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (71, 61)  => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (8, 9)    => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (3, 4)    => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  (5, 6)    => ([1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, …
  ⋮         => ⋮

This may seem strange, but we note that get_boundary_nodes(tri, pos), when pos is a vector of vertices, just returns pos once again so that get_boundary_nodes(get_boundary_nodes(tri, pos), ℓ) == u is always true, regardless of the boundary type. For example:

pos, ℓ = get_boundary_edge_map(tri3, 1, 2)
bn = get_boundary_nodes(tri3, pos)
37-element Vector{Int64}:
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
  ⋮
 81
 71
 61
 51
 41
 31
 21
 11
  1

For the sectioned case, we have:

tri4 = triangulate_rectangle(0, 1, 0, 1, 10, 10, single_boundary = false)
get_boundary_edge_map(tri4)
Dict{Tuple{Int64, Int64}, Tuple{Int64, Int64}} with 36 entries:
  (91, 81)  => (4, 1)
  (4, 5)    => (1, 4)
  (1, 2)    => (1, 1)
  (96, 95)  => (3, 5)
  (60, 70)  => (2, 6)
  (9, 10)   => (1, 9)
  (98, 97)  => (3, 3)
  (21, 11)  => (4, 8)
  (100, 99) => (3, 1)
  (95, 94)  => (3, 6)
  (93, 92)  => (3, 8)
  (30, 40)  => (2, 3)
  (31, 21)  => (4, 7)
  (61, 51)  => (4, 4)
  (7, 8)    => (1, 7)
  (71, 61)  => (4, 3)
  (8, 9)    => (1, 8)
  (3, 4)    => (1, 3)
  (5, 6)    => (1, 5)
  ⋮         => ⋮

When there is no constrained boundary, the boundary edge map will be empty:

tri5 = triangulate(rand(2, 50))
get_boundary_edge_map(tri5)
Dict{Tuple{Int64, Int64}, Tuple{Vector{Int64}, Int64}}()

get_ghost_vertex_map(tri)

This field is used to identify to what curve and to what section a ghost vertex is associated. For example, we have

get_ghost_vertex_map(tri)
Dict{Int64, Tuple{Int64, Int64}} with 8 entries:
  -5 => (4, 1)
  -1 => (1, 1)
  -8 => (6, 1)
  -7 => (5, 1)
  -3 => (2, 1)
  -2 => (1, 2)
  -4 => (3, 1)
  -6 => (4, 2)

We see, for example, that the ghost vertex -7 comes from the first section of the fifth boundary curve. In general, the mappings are of the form g => pos, where g is the ghost vertex and pos is such that get_boundary_nodes(tri, pos) is the boundary curve associated with g. For example, using map_ghost_vertex,

pos = map_ghost_vertex(tri, -7)
bn = get_boundary_nodes(tri, pos)
33-element Vector{Int64}:
   12
 1952
  415
  174
 1708
  364
   65
  616
  305
  544
    ⋮
   10
  551
  196
  586
   66
  611
   67
  558
   12

The forms of pos for the case of a contiguous boundary, a sectioned boundary, and no boundary are the same as for get_boundary_edge_map. For example,

get_ghost_vertex_map(tri3) # contiguous
Dict{Int64, Vector{Int64}} with 1 entry:
  -1 => [1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  91, 81, 71, 61, 51, 41, 31, 21, 11, …
get_ghost_vertex_map(tri4) # sectioned
Dict{Int64, Int64} with 4 entries:
  -1 => 1
  -3 => 3
  -2 => 2
  -4 => 4
get_ghost_vertex_map(tri5) # no boundary
Dict{Int64, Vector{Int64}} with 1 entry:
  -1 => []

get_ghost_vertex_ranges(tri)

This field is used to identify all the ghost vertices associated with a curve that has a specific ghost vertex. For example, we have

get_ghost_vertex_ranges(tri)
Dict{Int64, UnitRange{Int64}} with 8 entries:
  -5 => -6:-5
  -1 => -2:-1
  -8 => -8:-8
  -7 => -7:-7
  -3 => -3:-3
  -2 => -2:-1
  -4 => -4:-4
  -6 => -6:-5

This output means, for example, that the ghost vertex -5 is associated with a curve that has both ghost vertices -6 and -5 associated with it. You would use get_ghost_vertex_range to get the range of ghost vertices associated with a specific ghost vertex. For example,

get_ghost_vertex_range(tri, -5)
-6:-5

Similarly, the output for the contiguous, sectioned, and no boundary cases are

get_ghost_vertex_ranges(tri3) # contiguous
Dict{Int64, UnitRange{Int64}} with 1 entry:
  -1 => -1:-1
get_ghost_vertex_ranges(tri4) # sectioned
Dict{Int64, UnitRange{Int64}} with 4 entries:
  -1 => -4:-1
  -3 => -4:-1
  -2 => -4:-1
  -4 => -4:-1
get_ghost_vertex_ranges(tri5) # no boundary
Dict{Int64, UnitRange{Int64}} with 1 entry:
  -1 => -1:-1

Other Fields

There are some other more general fields to inspect.

get_convex_hull(tri)

For all triangulations, the convex_hull field stores the ConvexHull of the triangulation. For our triangulation, we have

get_convex_hull(tri)
Convex hull.
   Vertices:
24-element Vector{Int64}:
  84
  16
 105
  62
  51
  50
  49
  18
   1
   7
   ⋮
  53
  52
   6
   3
  17
  48
  60
  83
  84

To inspect the vertices of the convex hull, you use get_convex_hull_vertices:

get_convex_hull_vertices(tri)
24-element Vector{Int64}:
  84
  16
 105
  62
  51
  50
  49
  18
   1
   7
   ⋮
  53
  52
   6
   3
  17
  48
  60
  83
  84

If you ever need to reconstruct the convex hull, say after some dynamic updates, you would use convex_hull!.

get_representative_point_list(tri)

This field is related to the need for a representative point for each curve, as described in the ghost vertices section. For our triangulation, we have

DelaunayTriangulation.get_representative_point_list(tri)
Dict{Int64, DelaunayTriangulation.RepresentativeCoordinates{Int64, Float64}} with 6 entries:
  5 => RepresentativeCoordinates{Int64, Float64}(0.2, -2.9, 0)
  4 => RepresentativeCoordinates{Int64, Float64}(0.0654276, -3.01177, 0)
  6 => RepresentativeCoordinates{Int64, Float64}(0.0, -3.0, 0)
  2 => RepresentativeCoordinates{Int64, Float64}(8.65974e-15, 0.275, 0)
  3 => RepresentativeCoordinates{Int64, Float64}(0.0, -2.55, 0)
  1 => RepresentativeCoordinates{Int64, Float64}(-0.75, 0.25, 0)

This shows, for example, that the sixth boundary curve's representative point, i.e. its pole of inaccessibility, is the point (0, -3). You could also review this using DelaunayTriangulation.get_representative_point_coordinates:

DelaunayTriangulation.get_representative_point_coordinates(tri, 6)
(0.0, -3.0)

In the case of a triangulation with no constrained boundary, the representative point list is simply the centroid of the domain. If you ever need to recompute the representative points, you need DelaunayTriangulation.compute_representative_points!.

get_polygon_hierarchy(tri)

This field is used to store information about which boundary curves are contained in other boundary curves. For our triangulation, we have

DelaunayTriangulation.get_polygon_hierarchy(tri)
DelaunayTriangulation.PolygonHierarchy{Int64}(Bool[1, 0, 1, 1, 0, 0], DelaunayTriangulation.BoundingBox[[-2.0000000596046448, 2.0000000596046448] × [-7.450580596923828e-9, 0.5000000074505806], [-0.577350284704563, 0.5773502847045628] × [0.14999999627470972, 0.4000000037252903], [-2.0000000596046448, 2.0000000596046448] × [-5.0000000730156895, -0.09999992698431015], [-1.0000000298023222, 1.0000000124171444] × [-4.000000008207955, -1.9999999701976778], [-0.5000000208616256, 0.9000000208616257] × [-3.5000000178813933, -2.2999999821186066], [-1.1000000327825548, 1.1000000327825548] × [-4.100000032782554, -1.8999999672174452]], Dict{Int64, DelaunayTriangulation.PolygonTree{Int64}}(3 => PolygonTree at height 0 with index 3 and 1 children, 1 => PolygonTree at height 0 with index 1 and 1 children), DelaunayTriangulation.PolygonTree{Int64}[PolygonTree at height 2 with index 4 and 1 children])

This field is not intended for public use. One useful method that makes direct use of the polygon hierarchy is find_polygon, but the hierarchy's main use is in boundary enrichment for curve-bounded domains, and for checking arguments via check_args.

boundary_enricher(tri)

This field, just as for the polygon_hierarchy, is not intended for public use. It is the field used for enriching the boundary for initialising the triangulation of a curve-bounded domain.

DelaunayTriangulation.get_boundary_enricher(tri)
BoundaryEnricher with 8 boundary curves and 2068 points

get_cache(tri)

This field is not intended for public use. It is used to provide several caches for use during triangulation, such as reusing arrays for constrained triangulations.

DelaunayTriangulation.get_cache(tri)
TriangulationCache with storage.