DelaunayTriangulation.jl

Delaunay triangulations and Voronoi tessellations in two dimensions.
Author JuliaGeometry
Popularity
62 Stars
Updated Last
3 Months Ago
Started In
September 2022

DelaunayTriangulation

Stable Dev Coverage DOI Aqua QA

This is a package for constructing Delaunay triangulations and Voronoi tessellations of planar point sets. Supports unconstrained and constrained triangulations, weighted triangulations, mesh refinement, triangulation of curve bounded domains, Voronoi tessellations, power diagrams, and clipped and centroidal Voronoi tessellations. To install the package, do

julia>] add DelaunayTriangulation

Many features are available, some of these being:

  • Unconstrained and constrained Delaunay triangulations, supporting many types of domains.
  • Weighted Delaunay triangulations.
  • Computation of Voronoi tessellations, clipped Voronoi tessellations where the Voronoi tiles get clipped to a convex polygon, and centroidal Voronoi tessellations where each Voronoi tile's generator is the tile's centroid.
  • Power diagrams.
  • Mesh refinement, with support for custom angle and area constraints, as well as refinement of curve-bounded domains.
  • Dynamic point insertion, point deletion, and segment insertion, amongst many other operations.
  • Computation of convex hulls.
  • Triangulation of convex polygons.
  • Point location.
  • Computation of the pole of inaccessibility.
  • The interface for defining geometric primitives is fully customisable.

To ensure that the algorithms are robust, by default we make use of AdaptivePredicates.jl to use adaptive arithmetic for all geometric predicates in this package. This choice can be configured, allowing for the additional choices of ExactPredicates.jl or no adaptive or exact arithmetic at all; see the documentation. Much of the work in this package is derived from the book Delaunay Mesh Generation by Cheng, Dey, and Shewchuk (2013). Please see the documentation for much more information.

Some examples are below (and in the documentation), but if you would also like to see how DelaunayTriangulation.jl is used in other packages, see FiniteVolumeMethod.jl (solving 2D PDEs) and NaturalNeighbours.jl (scattered data interpolation).

If you would like to make an issue or contribute, please see the CONTRIBUTING.md file. For feature requests, please see the discussions tab.

Example Usage

Here is a quick example of some ways the package can be used. As mentioned, see the docs for many examples.

# using Pkg; Pkg.add(["DelaunayTriangulation", "CairoMakie", "StableRNGs"])
using DelaunayTriangulation, CairoMakie, StableRNGs

# Unconstrained 
points = rand(2, 50)
tri1 = triangulate(points) # default predicate kernel is AdaptiveKernel()

# Voronoi example 
vorn2 = voronoi(tri1)

# Clipped Voronoi 
vorn3 = voronoi(tri1, clip=true, predicates = ExactKernel())

# Smoothed Voronoi 
vorn4 = centroidal_smooth(vorn3; predicates = FastKernel()) # or do voronoi(tri1, clip = true, smooth = true)

# Constrained example with refinement 
boundary_points = [(0.0, 0.0), (1.0, 0.0), (1.0, 0.3), (0.5, 0.3),
    (0.3, 0.7), (0.1, 1.0), (0.0, 1.0), (0.0, 0.0)]
boundary_nodes, points = convert_boundary_points_to_indices(boundary_points)
tri5 = triangulate(points; boundary_nodes)
refine!(tri5; max_area=1e-2get_area(tri5))

# Disjoint constrained example with refinement 
boundary_points = [
    [[(0.0, 0.0), (1.0, 0.0), (1.0, 1.0), (0.0, 1.0), (0.0, 0.0)]],
    [[(0.3, 0.3), (0.3, 0.7), (0.7, 0.7), (0.7, 0.3), (0.3, 0.3)]],
    [[(1.2, 0.0), (1.4, 0.0), (1.4, 1.2), (0.0, 1.2), (0.0, 1.1),
        (1.1, 1.1), (1.1, 0.0), (1.2, 0.0)]]
]
boundary_nodes, points = convert_boundary_points_to_indices(boundary_points)
tri6 = triangulate(points; boundary_nodes)
refine!(tri6; max_area=1e-2get_area(tri6))

# Curve-bounded example
using DelaunayTriangulation: EllipticalArc
ellipse = EllipticalArc((1.0, 0.0), (1.0, 0.0), (0.0, 0.0), 1.0, 2.0, 0.0)
tri7 = triangulate(NTuple{2,Float64}[]; boundary_nodes=[ellipse])
refine!(tri7; max_area=1e-2get_area(tri7))

# Disjoint curve-bounded example 
ellipse = EllipticalArc((7.0, 3.5), (7.0, 3.5), (0.0, 3.5), 7.0, 10.0, 0.0)
catmull_cp = [(0.0, 0.0), (-2.0, -1.0), (-4.0, 0.0), (-5.0, 2.0), (-1.0, 4.0), (0.0, 3.0),
    (1.0, 4.0), (5.0, 2.0), (4.0, 0.0), (2.0, -1.0), (0.0, 0.0)]
catmull_spl = CatmullRomSpline(catmull_cp)
circle = CircularArc((0.5, 1.5), (0.5, 1.5), (0.0, 1.0))
circle2 = CircularArc((0.1, 1.5), (0.1, 1.5), (0.0, 1.0), positive=false)
points = [(-4.0, -10.0), (4.0, -10.0), (4.0, -7.0), (-4.0, -7.0)]
square = [1, 2, 3, 4, 1]
boundary_nodes = [[square], [[ellipse]], [[catmull_spl]], [[circle]], [[circle2]]]
tri8 = triangulate(points; boundary_nodes)
refine!(tri8; max_area=1e-2get_area(tri8)) # could also use find_polygon to help define a custom refinement function for each shape

# Weighted triangulation example
rng = StableRNG(1) # You can pass rngs for reproducibility. Alternatively, just use Random.seed!(n)
points = tuple.(rand(rng, 20), rand(rng, 20))
weights = 3randn(rng, 20)
tri9 = triangulate(points; weights, rng)

# Power diagram example 
points = tuple.(rand(rng, 100), rand(rng, 100))
weights = rand(rng, 100)
tri10 = triangulate(points; weights, rng)
vorn10 = voronoi(tri10) # can also use clip/smooth here 

# Clipped Voronoi example with a generic convex polygon
points = 10randn(rng, 2, 100)
weights = rand(rng, 100)
circ = CircularArc((0.0, 5.0), (0.0, 5.0), (0.0, 0.0)) # clip to a circle 
clip_points = [circ(t) for t in LinRange(0, 1, 100)]
clip_vertices = [1:(length(clip_points)-1); 1]
tri11 = triangulate(points; weights, rng)
vorn11 = voronoi(tri11, clip=true, clip_polygon=(clip_points, clip_vertices), rng=rng)

# Plotting 
fig = Figure(fontsize = 42, size = (2800, 2200))
wh = (width = 600, height = 600)
ax = Axis(fig[1, 1]; title = "Unconstrained", wh...);            triplot!(ax, tri1)
ax = Axis(fig[1, 2]; title = "Voronoi", wh...);                  voronoiplot!(ax, vorn2)
ax = Axis(fig[1, 3]; title = "Clipped Voronoi", wh...);          voronoiplot!(ax, vorn3)
ax = Axis(fig[1, 4]; title = "Centroidal Voronoi", wh...);       voronoiplot!(ax, vorn4)
ax = Axis(fig[2, 1]; title = "Constrained", wh...);              triplot!(ax, tri5)
ax = Axis(fig[2, 2]; title = "Disjoint Constrained", wh...);     triplot!(ax, tri6)
ax = Axis(fig[2, 3]; title = "Curve-Bounded", wh...);            triplot!(ax, tri7)
ax = Axis(fig[2, 4]; title = "Disjoint Curve-Bounded", wh...);   triplot!(ax, tri8)
ax = Axis(fig[3, 1]; title = "Weighted", wh...);                 triplot!(ax, tri9)
ax = Axis(fig[3, 2]; title = "Power Diagram", wh...);            voronoiplot!(ax, vorn10)
ax = Axis(fig[3, 3]; title = "Generic Clipped Voronoi", wh...);  voronoiplot!(ax, vorn11)

Similar Packages

This is not the only Delaunay triangulation package available. Some others are:

  • VoronoiDelaunay.jl: A pure Julia library that constructs planar triangulations and tessellations like in this package, although no support for constrained triangulations / mesh refinement or clipped / centroid tessellations. Restricts points to $[1, 2] \times [1, 2]$.
  • VoronoiCells.jl: A pure Julia library that extends VoronoiDelaunay.jl. This package provides useful tools for constructing and working with Voronoi tessellations. Supports clipping Voronoi cells to a specified rectangle. Like VoronoiDelaunay.jl, restricts points to $[1, 2] \times [1, 2]$.
  • Delaunay.jl: Wraps Python's main Delaunay triangulation library, scipy.spatial.Delaunay, for computing Delaunay triangulations in $\mathbb R^N$. I don't believe constrained triangulations or mesh refinement is available here.
  • MiniQhull.jl: Wraps Qhull for computing unconstrained Delaunay triangulations in $\mathbb R^N$. No support is provided for mesh refinement.
  • DirectQhull.jl: Similar to MiniQhull.jl, although also provides support for convex hulls and Voronoi tessellations from Qhull.
  • Delaunator.jl: A pure Julia library modelled after the JavaScript Delaunator library. This package can construct unconstrained triangulations of planar point sets. No support is available for constrained triangulations or mesh refinement, although support exists for computing the dual Voronoi tessellation. Centroidal tessellations are not implemented, although the Voronoi cells can be clipped to a bounding box.
  • TriangleMesh.jl, Triangulate.jl, Triangle.jl: Interfaces to Shewchuk's Triangle library.
  • TetGen.jl: This is for Delaunay tetrahedralisation, wrapping TetGen.
  • GMT.jl: A wrapper of GMT, allowing for unconstrained Delaunay triangulations in two dimensions, and for spherical triangulation, i.e. triangulation of points lying on a sphere.
  • Quickhull.jl: A pure Julia library for unconstrained triangulations, Voronoi tessellations, and convex hulls in $N$ dimensions.