pub struct DelaunayTriangulation<V, DE = (), UE = (), F = (), L = LastUsedVertexHintGenerator>
where V: HasPosition, DE: Default, UE: Default, F: Default, L: HintGenerator<<V as HasPosition>::Scalar>,
{ /* private fields */ }
Expand description

A two dimensional Delaunay triangulation.

A Delaunay triangulation a triangulation that fulfills the Delaunay Property: No vertex of the triangulation is contained in the circumcircle of any triangle. As a consequence, Delaunay triangulations are well suited to support interpolation algorithms and nearest neighbor searches. It is often constructed in order to retrieve its dual graph, the Voronoi diagram.

An example triangulation with a few circumcircles drawn. No circumcircle contains more than three vertices.

Most methods on this type require the Triangulation trait. Refer to its documentation for more details on how to use DelaunayTriangulation.

Basic Usage

Vertices need to implement the HasPosition trait. Spade bundles the Point2 struct for basic use cases.

Basic example

use spade::{DelaunayTriangulation, Triangulation, Point2, InsertionError};

fn main() -> Result<(), InsertionError> {

   let mut triangulation: DelaunayTriangulation<_> = DelaunayTriangulation::new();

   // Insert three vertices that span one triangle (face)
   triangulation.insert(Point2::new(0.0, 1.0))?;
   triangulation.insert(Point2::new(1.0, 1.0))?;
   triangulation.insert(Point2::new(0.5, -1.0))?;

   assert_eq!(triangulation.num_vertices(), 3);
   assert_eq!(triangulation.num_inner_faces(), 1);
   assert_eq!(triangulation.num_undirected_edges(), 3);
   Ok(())
}

Right handed and left handed coordinate systems

For simplicity, all method names and their documentation assume that the underlying coordinate system is right handed (e.g. x axis points to the right, y axis points upwards). If a left handed system (lhs) is used, any term related to orientation needs to be reversed:

  • “left” becomes “right” (example: the face of a directed edge is on the right side for a lhs
  • “counter clock wise” becomes “clockwise” (example: the vertices of a face are returned in clock wise order for a lhs)
left handed systemright handed system
# Extracting geometry information

Spade uses handles to extract the triangulation’s geometry. Handles are usually retrieved by inserting a vertex or by iterating.

Example

 fn main() -> Result<(), spade::InsertionError> {
use crate::spade::{DelaunayTriangulation, Triangulation, Point2};

let mut triangulation: DelaunayTriangulation<Point2<f64>> = DelaunayTriangulation::new();

triangulation.insert(Point2::new(0.0, 1.0))?;
triangulation.insert(Point2::new(1.0, 1.0))?;
triangulation.insert(Point2::new(0.5, -1.0))?;

for face in triangulation.inner_faces() {
  // face is a FaceHandle
  // edges is an array containing 3 directed edge handles
  let edges = face.adjacent_edges();
  for edge in &edges {
    let from = edge.from();
let to = edge.to();
    // from and to are vertex handles
    println!("found an edge: {:?} -> {:?}", from, to);
  }

  // vertices is an array containing 3 vertex handles
  let vertices = face.vertices();
  for vertex in &vertices {
    println!("Found vertex with position {:?}", vertex.position());
  }
}

Type parameters

The triangulation’s vertices, edges and faces can contain custom data. By default, the edge and face types are set to (). The vertex type must be specified.

  • V: HasPosition The vertex type
  • DE: Default The directed edge type.
  • UE: Default The undirected edge type.
  • F: Default The face type.

Only vertices can be inserted directly. Faces and edges are create via Default::default(). Usually, edge and face data will need to be modified in a separate pass.

Setting any custom data works by calling vertex_data_mut, directed_edge_data_mut, undirected_edge_data_mut and face_data_mut.

Example

fn main() -> Result<(), spade::InsertionError> {
use crate::spade::{DelaunayTriangulation, Triangulation, Point2};

// A custom undirected edge type used to cache the length of an edge
#[derive(Default)]
struct EdgeWithLength { length: f64 }

// Creates a new triangulation with a custom undirected edge type
let mut triangulation: DelaunayTriangulation<Point2<f64>, (), EdgeWithLength>
                         = DelaunayTriangulation::new();

triangulation.insert(Point2::new(0.0, 1.0))?;
triangulation.insert(Point2::new(1.0, 1.0))?;
triangulation.insert(Point2::new(0.5, -1.0))?;

for edge in triangulation.fixed_undirected_edges() {
  let positions = triangulation.undirected_edge(edge).positions();
  let length = positions[0].distance_2(positions[1]).sqrt();
  // Write length into the edge data
  triangulation.undirected_edge_data_mut(edge).length = length;
}

for edge in triangulation.undirected_edges() {
   let length = edge.data().length;
   assert!(length > 0.0);
}

Outer face

Every triangulation contains an outer face that is adjacent to all edges of the triangulation’s convex hull. The outer face is even present for a triangulation without vertices. It is either referenced by calling Triangulation::outer_face() or handles::OUTER_FACE

inner outer inner outer inner outer

Voronoi Diagram

the dual graph of the Delaunay triangulation is the Voronoi diagram. The Voronoi diagram subdivides the plane into several Voronoi cells (called VoronoiFace by Spade for consistency) which contain all points in the plane that share a common nearest neighbor.

Keep in mind that “faces” and “vertices” are swapped - an (inner) Voronoi vertex corresponds to a single Delaunay face. The position of an inner voronoi vertex is the circumcenter of its dual Delaunay face.

A Delaunay triangulation (black lines) and its dual graph, the Voronoi diagram (blue lines)

Extracting the Voronoi Diagram

Spade defines various functions to extract information about the Voronoi diagram:

Types

Iterators

Conversion

Extracting the Voronoi Diagram (Example)

Extracting the geometry of the voronoi diagram can be slightly tricky as some of the voronoi extend into infinity (see the dashed lines in the example above).

use spade::handles::{VoronoiVertex::Inner, VoronoiVertex::Outer};
use spade::{DelaunayTriangulation, Point2, Triangulation};

// Prints out the location of all voronoi edges in a triangulation
fn log_voronoi_diagram(triangulation: &DelaunayTriangulation<Point2<f64>>) {
    for edge in triangulation.undirected_voronoi_edges() {
        match edge.vertices() {
            [Inner(from), Inner(to)] => {
                // "from" and "to" are inner faces of the Delaunay triangulation
                println!(
                    "Found voronoi edge between {:?} and {:?}",
                    from.circumcenter(),
                    to.circumcenter()
                );
            }
            [Inner(from), Outer(edge)] | [Outer(edge), Inner(from)] => {
                // Some lines don't have a finite end and extend into infinity.
                println!(
                    "Found infinite voronoi edge going out of {:?} into the direction {:?}",
                    from.circumcenter(),
                    edge.direction_vector()
                );
            }
            [Outer(_), Outer(_)] => {
                // This case only happens if all vertices of the triangulation lie on the
                // same line and can probably be ignored.
            }
        }
    }
}

Performance tuning

Fine-tuning a Delaunay triangulation can be more tricky from time to time. However, some will nearly always be the right thing to do:

  • Measure, don’t guess. Execution times are hard to predict.
  • If you plan to perform several random access queries (e.g. looking up the point at an arbitrary position): Consider using `HierarchyHintGenerator
  • For data sets with uniformly distributed vertices: Use HierarchyHintGenerator if bulk loading is not applicable.
  • For data sets where vertices are inserted in close local proximity (each vertex is not too far away from the previously inserted vertex): Consider using LastUsedVertexHintGenerator.
  • Try to avoid large custom data types for edges, vertices and faces.
  • Using f64 and f32 as scalar type will usually end up roughly having the same run time performance.
  • Prefer using bulk_load over insert.
  • The run time of all vertex operations (insertion, removal and lookup) is roughly the same for larger triangulations.

Complexity classes

This table display the average and amortized cost for inserting a vertex into a triangulation with n vertices.

Uniformly distributed verticesInsertion of vertices with local proximity
LastUsedVertexHintGeneratorO(sqrt(n)) (worst case)O(1) (best case), fastest
HierarchyHintGeneratorO(log(n)) (Average case)O(1) (best case)

See also

Delaunay triangulations are closely related to constrained Delaunay triangulations

Implementations§

source§

impl<V, DE, UE, F, L> DelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition, DE: Default, UE: Default, F: Default, L: HintGenerator<<V as HasPosition>::Scalar>,

source

pub fn nearest_neighbor( &self, position: Point2<<V as HasPosition>::Scalar> ) -> Option<VertexHandle<'_, V, DE, UE, F>>

Returns the nearest neighbor of a given input vertex.

Returns None if the triangulation is empty.

Runtime

This method takes O(sqrt(n)) on average where n is the number of vertices.

source

pub fn bulk_load_stable(elements: Vec<V>) -> Result<Self, InsertionError>

Creates a new delaunay triangulation with an efficient bulk loading strategy.

In contrast to Triangulation::bulk_load, this method will create a triangulation with vertices returned in the same order as the input vertices.

Duplicate handling

If two vertices have the same position, only one of them will be included in the final triangulation. It is undefined which of them is discarded.

For example, if the input vertices are [v0, v1, v2, v1] (where v1 is duplicated), the resulting triangulation will be either [v0, v1, v2] or [v0, v2, v1]

Consider checking the triangulation’s size after calling this method to ensure that no duplicates were present.

Example
use spade::{DelaunayTriangulation, Point2, Triangulation};

let vertices = vec![
     Point2::new(0.5, 1.0),
     Point2::new(-0.5, 2.0),
     Point2::new(0.2, 3.0),
     Point2::new(0.0, 4.0),
     Point2::new(-0.2, 5.0)
];
let triangulation = DelaunayTriangulation::<Point2<f64>>::bulk_load_stable(vertices.clone())?;
// This assert will not hold for regular bulk loading!
assert_eq!(triangulation.vertices().map(|v| *v.data()).collect::<Vec<_>>(), vertices);

// This is how you would check for duplicates:
let duplicates_found = triangulation.num_vertices() < vertices.len();
assert!(!duplicates_found);
source§

impl<V, DE, UE, F, L> DelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition, DE: Default, UE: Default, F: Default, V::Scalar: Float, L: HintGenerator<<V as HasPosition>::Scalar>,

source

pub fn natural_neighbor(&self) -> NaturalNeighbor<'_, Self>

Allows using natural neighbor interpolation on this triangulation. Refer to the documentation of NaturalNeighbor for more information.

Trait Implementations§

source§

impl<V, DE, UE, F, L> Clone for DelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition + Clone, DE: Default + Clone, UE: Default + Clone, F: Default + Clone, L: HintGenerator<<V as HasPosition>::Scalar> + Clone,

source§

fn clone(&self) -> DelaunayTriangulation<V, DE, UE, F, L>

Returns a copy of the value. Read more
1.0.0 · source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
source§

impl<V, DE, UE, F, L> Debug for DelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition + Debug, DE: Default + Debug, UE: Default + Debug, F: Default + Debug, L: HintGenerator<<V as HasPosition>::Scalar> + Debug,

source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
source§

impl<V, DE, UE, F, L> Default for DelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition, DE: Default, UE: Default, F: Default, L: HintGenerator<<V as HasPosition>::Scalar>,

source§

fn default() -> Self

Returns the “default value” for a type. Read more
source§

impl<V, DE, UE, F, L> From<DelaunayTriangulation<V, DE, UE, F, L>> for ConstrainedDelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition, DE: Default, UE: Default, F: Default, L: HintGenerator<<V as HasPosition>::Scalar>,

source§

fn from(value: DelaunayTriangulation<V, DE, UE, F, L>) -> Self

Converts to this type from the input type.
source§

impl<V, DE, UE, F, L> Triangulation for DelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition, DE: Default, UE: Default, F: Default, L: HintGenerator<<V as HasPosition>::Scalar>,

§

type Vertex = V

The triangulation’s vertex type.
§

type DirectedEdge = DE

The triangulation’s edge type. Any new edge is created by using the Default trait.
§

type UndirectedEdge = UE

The triangulation’s undirected edge type. Any new edge is created by using the Default trait.
§

type Face = F

The triangulation’s face type. Any new face is created by using the Default trait.
§

type HintGenerator = L

The hint generator used by the triangulation. See HintGenerator for more information.
source§

fn new() -> Self

Creates a new triangulation. Read more
source§

fn with_capacity( num_vertices: usize, num_undirected_edges: usize, num_faces: usize ) -> Self

Creates a new triangulation and pre-allocates some space for vertices, edges and faces
source§

fn clear(&mut self)

Removes all edges, faces and vertices from the triangulation. Read more
source§

fn bulk_load(elements: Vec<Self::Vertex>) -> Result<Self, InsertionError>

Creates a new triangulation populated with some vertices. Read more
source§

fn vertex( &self, handle: FixedVertexHandle ) -> VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

Converts a fixed vertex handle to a reference vertex handle. Read more
source§

fn vertex_data_mut(&mut self, handle: FixedVertexHandle) -> &mut Self::Vertex

Returns a mutable reference to the associated data of a vertex.
source§

fn face<InnerOuter: InnerOuterMarker>( &self, handle: FixedFaceHandle<InnerOuter> ) -> FaceHandle<'_, InnerOuter, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

Converts a fixed face handle to a reference face handle. Read more
source§

fn outer_face( &self ) -> FaceHandle<'_, PossiblyOuterTag, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

Returns a reference handle to the single outer face of this triangulation.
source§

fn directed_edge( &self, handle: FixedDirectedEdgeHandle ) -> DirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

Converts a fixed directed edge handle handle to a reference directed edge handle. Read more
source§

fn undirected_edge( &self, handle: FixedUndirectedEdgeHandle ) -> UndirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

Converts a fixed undirected edge handle to a reference undirected edge handle. Read more
source§

fn undirected_edge_data_mut( &mut self, handle: FixedUndirectedEdgeHandle ) -> &mut Self::UndirectedEdge

Returns a mutable reference ot the associated data of an undirected edge.
source§

fn num_all_faces(&self) -> usize

Returns the number of all faces, including the single outer face, of this triangulation. Read more
source§

fn num_inner_faces(&self) -> usize

Returns the number of inner faces in this triangulation.
source§

fn num_undirected_edges(&self) -> usize

Returns the number of undirected edges in this triangulation.
source§

fn num_directed_edges(&self) -> usize

Returns the number of directed edges in this triangulation.
source§

fn convex_hull_size(&self) -> usize

Returns the number of edges of the convex hull. Read more
source§

fn directed_edges( &self ) -> DirectedEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all directed edges. Read more
source§

fn undirected_edges( &self ) -> UndirectedEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator over all undirected edges. Read more
source§

fn num_vertices(&self) -> usize

Returns the number vertices in this triangulation.
source§

fn vertices( &self ) -> VertexIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all vertices. Read more
source§

fn fixed_vertices(&self) -> FixedVertexIterator

An iterator visiting all vertices. Read more
source§

fn all_faces( &self ) -> FaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all faces. Read more
source§

fn inner_faces( &self ) -> InnerFaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all inner faces. The iterator type is FaceHandle<InnerTag, …>.
source§

fn voronoi_faces( &self ) -> VoronoiFaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all faces of the Voronoi diagram. Read more
source§

fn directed_voronoi_edges( &self ) -> DirectedVoronoiEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all directed voronoi edges. Read more
source§

fn undirected_voronoi_edges( &self ) -> UndirectedVoronoiEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

An iterator visiting all undirected voronoi edges. Read more
source§

fn locate( &self, point: Point2<<Self::Vertex as HasPosition>::Scalar> ) -> PositionInTriangulation

Returns information about the location of a point in a triangulation.
source§

fn locate_vertex( &self, point: Point2<<Self::Vertex as HasPosition>::Scalar> ) -> Option<VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>

Locates a vertex at a given position. Read more
source§

fn get_edge_from_neighbors( &self, from: FixedVertexHandle, to: FixedVertexHandle ) -> Option<DirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>

Returns an edge between two vertices. Read more
source§

fn locate_with_hint( &self, point: Point2<<Self::Vertex as HasPosition>::Scalar>, hint: FixedVertexHandle ) -> PositionInTriangulation

Returns information about the location of a point in a triangulation. Read more
source§

fn insert_with_hint( &mut self, t: Self::Vertex, hint: FixedVertexHandle ) -> Result<FixedVertexHandle, InsertionError>

Inserts a new vertex into the triangulation. Read more
source§

fn locate_and_remove( &mut self, point: Point2<<Self::Vertex as HasPosition>::Scalar> ) -> Option<Self::Vertex>

Attempts to remove a vertex from the triangulation. Read more
source§

fn remove(&mut self, vertex: FixedVertexHandle) -> Self::Vertex

Removes a vertex from the triangulation. Read more
source§

fn insert( &mut self, vertex: Self::Vertex ) -> Result<FixedVertexHandle, InsertionError>

Inserts a new vertex into the triangulation. Read more
source§

fn fixed_undirected_edges(&self) -> FixedUndirectedEdgeIterator

An iterator visiting all undirected edges. Read more
source§

fn fixed_directed_edges(&self) -> FixedDirectedEdgeIterator

An iterator visiting all directed edges. Read more
source§

fn fixed_all_faces(&self) -> FixedFaceIterator

An iterator visiting all faces. Read more
source§

fn fixed_inner_faces(&self) -> FixedInnerFaceIterator

An iterator visiting all inner faces of the triangulation. Read more
source§

fn all_vertices_on_line(&self) -> bool

Returns true if all vertices lie on a single line. Read more
source§

fn convex_hull( &self ) -> HullIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>

Returns an iterator over all convex hull edges. Read more
source§

fn face_data_mut<InnerOuter: InnerOuterMarker>( &mut self, handle: FixedFaceHandle<InnerOuter> ) -> &mut Self::Face

Returns a mutable reference to the associated data of a face.
source§

fn directed_edge_data_mut( &mut self, handle: FixedDirectedEdgeHandle ) -> &mut Self::DirectedEdge

Returns a mutable reference to the associated data of a directed edge.

Auto Trait Implementations§

§

impl<V, DE, UE, F, L> RefUnwindSafe for DelaunayTriangulation<V, DE, UE, F, L>

§

impl<V, DE, UE, F, L> Send for DelaunayTriangulation<V, DE, UE, F, L>
where DE: Send, F: Send, L: Send, UE: Send, V: Send,

§

impl<V, DE, UE, F, L> Sync for DelaunayTriangulation<V, DE, UE, F, L>
where DE: Sync, F: Sync, L: Sync, UE: Sync, V: Sync,

§

impl<V, DE, UE, F, L> Unpin for DelaunayTriangulation<V, DE, UE, F, L>
where DE: Unpin, F: Unpin, L: Unpin, UE: Unpin, V: Unpin,

§

impl<V, DE, UE, F, L> UnwindSafe for DelaunayTriangulation<V, DE, UE, F, L>

Blanket Implementations§

source§

impl<T> Any for T
where T: 'static + ?Sized,

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

impl<T> Borrow<T> for T
where T: ?Sized,

source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
source§

impl<T> From<T> for T

source§

fn from(t: T) -> T

Returns the argument unchanged.

source§

impl<T, U> Into<U> for T
where U: From<T>,

source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

source§

impl<T> ToOwned for T
where T: Clone,

§

type Owned = T

The resulting type after obtaining ownership.
source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

§

type Error = Infallible

The type returned in the event of a conversion error.
source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.