src/delaunator

  Source   Edit

Delaunator-Nim is a port of Mapbox/Delaunator, a fast library for 2D Delaunay triangulation. In addition, this port includes a set of helpers for utilizing the structures of a Delaunator object, as well as clipping of infinite regions of the Voronoi diagram.

A Delaunator object is constructed from a set of points (aka 'sites'), and contains two key compact sequences of integers, triangles and halfedges. This representation of the triangulation, while less convenient, is what makes the library fast.

Delaunator construction can originate from a flat seq of coordinates,

Example:

import src/delaunator
var
  # A flat seq of `float64` coordinates
  coords = @[63.59410858154297, 198.1050262451172, 215.7989349365234, 171.0301208496094,
             33.8256950378418,  261.359130859375, 40.81229019165039, 61.88509368896484,
             189.6730651855469, 168.2080078125, 247.6787414550781, 222.6421508789062,
             265.9251403808594, 81.62255096435547, 21.60958862304688, 253.6200256347656,
             24.65586090087891, 67.60309600830078, 27.14787483215332, 113.4554977416992]
  # Construct
  d = delaunator.fromCoords[float64](coords)

# Triplets of site ids.
echo d.triangles
Output:
@[4, 5, 1, 4, 0, 5, 5, 6, 1, 1, 6, 4, 4, 9, 0, 0, 2, 5, 0, 7, 2, 3, 9, 4, 0, 9, 7, 6, 3, 4, 3, 8, 9, 9, 8, 7]

or from some pairwise sequence, with conversion to a float type,

Example:

import src/delaunator
var
  # A pairwise seq of `int` points
  points = @[
    [63, 198], [215, 171], [33,  261], [40, 61], [189, 168],
    [247, 222], [265, 81], [21, 253], [24, 67], [27, 113]
  ]
  # Construct into `float32` coordinates
  d = delaunator.fromPoints[array[2, int], float32](points)

# Halfedges of triangulation.
echo d.halfedges
Output:
@[5, 8, 11, 14, 17, 0, -1, 9, 1, 7, 29, 2, 22, 24, 3, 20, -1, 4, 26, -1, 15, 32, 12, 28, 13, 35, 18, -1, 23, 10, -1, 33, 21, 31, -1, 25]

Construction from custom types is aided by fromCustom.

Both fields, triangles and halfedges are sequences indexed by halfedge id. Importantly, notice that some halfedges index to '-1'. These halfedges have no opposite because they reside on the triangulation's hull. To quote Mapbox's guide to Delaunator's data structures:

A triangle edge may be shared with another triangle. Instead of thinking about each edge A ↔︎ B, we will use two half-edges A → B and B → A. Having two half-edges is the key to everything this library provides.

Half-edges e are the indices into both of delaunator’s outputs:

  • delaunay.triangles[e] returns the point id where the half-edge starts
  • delaunay.halfedges[e] returns the opposite half-edge in the adjacent triangle, or -1 if there is no adjacent triangle

Triangle ids and half-edge ids are related.

  • The half-edges of triangle t are 3 * t, 3 * t + 1, and 3 * t + 2.
  • The triangle of half-edge id e is floor(e/3).

The above linked guide is still very applicable to this port, and the helpers described therein, along with many more, are implemented in delaunator/helpers

Types

Delaunator[T] = ref object
  coords*: seq[T]            ## Flattened sequence of site points.
  minX*, minY*, maxX*, maxY*: T ## Extents of *coords*.
  triangles*: seq[uint32]    ## Sequence of triplet indices into *coords* defining delaunay triangulation.
  halfedges*: seq[int32]     ## Sequence of complement halfedges to that of the index.
  hull*: seq[uint32]         ## Sequence of point ids comprising the triangulation's hull.
  vectors*: seq[T]           ## Sequence of rays emanating from each triangle circumcenter adjacent to a hull site. Used for clipping infinite Voronoi regions.
  bounds*: tuple[minX, minY, maxX, maxY: T] ## Clipping bounds for the infinate Voronoi regions.
  trianglesLen: int32
  d_triangles: seq[uint32]
  d_halfedges: seq[int32]
  d_hashSize: int
  d_hullStart: int
  d_hullPrev: seq[uint32]
  d_hullNext: seq[uint32]
  d_hullTri: seq[uint32]
  d_hullHash: seq[int32]
  d_ids: seq[uint32]
  d_dists: seq[T]
  d_pointToLeftmostHalfedgeIndex: Table[uint32, int32]
This object holds the datastructures neccessary to build and navigate the Delaunay-Voronoi dual graph.   Source   Edit

Procs

func boundHeight[T](d: Delaunator[T]): T
Height of defined bounds.   Source   Edit
func boundWidth[T](d: Delaunator[T]): T
Width of defined bounds.   Source   Edit
func circumcenter[F](ax, ay, bx, by, cx, cy: F): tuple[x, y: F] {.inline.}
  Source   Edit
func extentHeight[T](d: Delaunator[T]): T
Height of the triangulation.   Source   Edit
func extentWidth[T](d: Delaunator[T]): T
Width of the triangulation.   Source   Edit
proc fromCoords[T](coordinates: var seq[T]): Delaunator[T]
Returns a Delaunator object constructed from coordinates, a flattened sequence of points representing site locations.   Source   Edit
proc fromCustom[P, T](points: seq[P]; getX, getY: proc (p: P): T): Delaunator[T]
Returns a Delaunator object constructed from points, a sequence of some type from which x and y coordinate values are extracted via the specified getX and getY procs.   Source   Edit
proc fromPoints[P, T](points: seq[P]): Delaunator[T]
Returns a Delaunator object constructed from points, a sequence of some pairwise type from which x and y coordinate values can be extracted via the '[]' operator. When x and y values are not at [0] and [1] respectively, use fromCustom instead.   Source   Edit
func hullNext(d: Delaunator; sid: uint32): uint32
Returns the id of the next site of the hull following the site defined by sid.   Source   Edit
func hullPrev(d: Delaunator; sid: uint32): uint32
Returns the id of the previous site of the hull preceding the site defined by sid.   Source   Edit
func siteToLeftmostHalfedge(d: Delaunator; sid: uint32): int32
Returns the id of the 'leftmost' incomming halfedge to the site defined by sid. 'Leftmost' can be understood as if one were standing at the site's position. Used for constructing Voronoi edges / regions.   Source   Edit
proc update[T](this: var Delaunator)
Builds and rebuilds the dual-graph. This procedure should be called on a Delaunator object after changing the object's coords.   Source   Edit