Struct spade::DelaunayTriangulation
source · 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 system | right handed system |
---|---|
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 typeDE: 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
§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
- DirectedVoronoiEdge::as_undirected()
- UndirectedVoronoiEdge::as_directed()
- DirectedEdgeHandle::as_voronoi_edge()
- DirectedVoronoiEdge::as_delaunay_edge()
- UndirectedEdgeHandle::as_voronoi_edge()
- UndirectedVoronoiEdge::as_delaunay_edge()
§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
andf32
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 vertices | Insertion of vertices with local proximity | |
---|---|---|
LastUsedVertexHintGenerator | O(sqrt(n)) (worst case) | O(1) (best case), fastest |
HierarchyHintGenerator | O(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>,
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>,
sourcepub fn nearest_neighbor(
&self,
position: Point2<<V as HasPosition>::Scalar>
) -> Option<VertexHandle<'_, V, DE, UE, F>>
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.
sourcepub fn bulk_load_stable(elements: Vec<V>) -> Result<Self, InsertionError>
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>,
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>,
sourcepub fn natural_neighbor(&self) -> NaturalNeighbor<'_, Self>
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,
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>
fn clone(&self) -> DelaunayTriangulation<V, DE, UE, F, L>
1.0.0 · source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source
. Read moresource§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,
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§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>,
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§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>,
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
fn from(value: DelaunayTriangulation<V, DE, UE, F, L>) -> Self
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>,
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 DirectedEdge = DE
type DirectedEdge = DE
Default
trait.§type UndirectedEdge = UE
type UndirectedEdge = UE
Default
trait.§type HintGenerator = L
type HintGenerator = L
source§fn with_capacity(
num_vertices: usize,
num_undirected_edges: usize,
num_faces: usize
) -> Self
fn with_capacity( num_vertices: usize, num_undirected_edges: usize, num_faces: usize ) -> Self
source§fn bulk_load(elements: Vec<Self::Vertex>) -> Result<Self, InsertionError>
fn bulk_load(elements: Vec<Self::Vertex>) -> Result<Self, InsertionError>
source§fn vertex(
&self,
handle: FixedVertexHandle
) -> VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn vertex( &self, handle: FixedVertexHandle ) -> VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn vertex_data_mut(&mut self, handle: FixedVertexHandle) -> &mut Self::Vertex
fn vertex_data_mut(&mut self, handle: FixedVertexHandle) -> &mut Self::Vertex
source§fn face<InnerOuter: InnerOuterMarker>(
&self,
handle: FixedFaceHandle<InnerOuter>
) -> FaceHandle<'_, InnerOuter, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn face<InnerOuter: InnerOuterMarker>( &self, handle: FixedFaceHandle<InnerOuter> ) -> FaceHandle<'_, InnerOuter, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn outer_face(
&self
) -> FaceHandle<'_, PossiblyOuterTag, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn outer_face( &self ) -> FaceHandle<'_, PossiblyOuterTag, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn directed_edge(
&self,
handle: FixedDirectedEdgeHandle
) -> DirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn directed_edge( &self, handle: FixedDirectedEdgeHandle ) -> DirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn undirected_edge(
&self,
handle: FixedUndirectedEdgeHandle
) -> UndirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn undirected_edge( &self, handle: FixedUndirectedEdgeHandle ) -> UndirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn undirected_edge_data_mut(
&mut self,
handle: FixedUndirectedEdgeHandle
) -> &mut Self::UndirectedEdge
fn undirected_edge_data_mut( &mut self, handle: FixedUndirectedEdgeHandle ) -> &mut Self::UndirectedEdge
source§fn num_all_faces(&self) -> usize
fn num_all_faces(&self) -> usize
source§fn num_inner_faces(&self) -> usize
fn num_inner_faces(&self) -> usize
source§fn num_undirected_edges(&self) -> usize
fn num_undirected_edges(&self) -> usize
source§fn num_directed_edges(&self) -> usize
fn num_directed_edges(&self) -> usize
source§fn convex_hull_size(&self) -> usize
fn convex_hull_size(&self) -> usize
source§fn directed_edges(
&self
) -> DirectedEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn directed_edges( &self ) -> DirectedEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn undirected_edges(
&self
) -> UndirectedEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn undirected_edges( &self ) -> UndirectedEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn num_vertices(&self) -> usize
fn num_vertices(&self) -> usize
source§fn vertices(
&self
) -> VertexIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn vertices( &self ) -> VertexIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn fixed_vertices(&self) -> FixedVertexIterator
fn fixed_vertices(&self) -> FixedVertexIterator
source§fn all_faces(
&self
) -> FaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn all_faces( &self ) -> FaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn inner_faces(
&self
) -> InnerFaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn inner_faces( &self ) -> InnerFaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn voronoi_faces(
&self
) -> VoronoiFaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn voronoi_faces( &self ) -> VoronoiFaceIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn directed_voronoi_edges(
&self
) -> DirectedVoronoiEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn directed_voronoi_edges( &self ) -> DirectedVoronoiEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn undirected_voronoi_edges(
&self
) -> UndirectedVoronoiEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
fn undirected_voronoi_edges( &self ) -> UndirectedVoronoiEdgeIterator<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>
source§fn locate(
&self,
point: Point2<<Self::Vertex as HasPosition>::Scalar>
) -> PositionInTriangulation
fn locate( &self, point: Point2<<Self::Vertex as HasPosition>::Scalar> ) -> PositionInTriangulation
source§fn locate_vertex(
&self,
point: Point2<<Self::Vertex as HasPosition>::Scalar>
) -> Option<VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>
fn locate_vertex( &self, point: Point2<<Self::Vertex as HasPosition>::Scalar> ) -> Option<VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>
source§fn get_edge_from_neighbors(
&self,
from: FixedVertexHandle,
to: FixedVertexHandle
) -> Option<DirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>
fn get_edge_from_neighbors( &self, from: FixedVertexHandle, to: FixedVertexHandle ) -> Option<DirectedEdgeHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>
source§fn locate_with_hint(
&self,
point: Point2<<Self::Vertex as HasPosition>::Scalar>,
hint: FixedVertexHandle
) -> PositionInTriangulation
fn locate_with_hint( &self, point: Point2<<Self::Vertex as HasPosition>::Scalar>, hint: FixedVertexHandle ) -> PositionInTriangulation
source§fn insert_with_hint(
&mut self,
t: Self::Vertex,
hint: FixedVertexHandle
) -> Result<FixedVertexHandle, InsertionError>
fn insert_with_hint( &mut self, t: Self::Vertex, hint: FixedVertexHandle ) -> Result<FixedVertexHandle, InsertionError>
source§fn locate_and_remove(
&mut self,
point: Point2<<Self::Vertex as HasPosition>::Scalar>
) -> Option<Self::Vertex>
fn locate_and_remove( &mut self, point: Point2<<Self::Vertex as HasPosition>::Scalar> ) -> Option<Self::Vertex>
source§fn remove(&mut self, vertex: FixedVertexHandle) -> Self::Vertex
fn remove(&mut self, vertex: FixedVertexHandle) -> Self::Vertex
source§fn insert(
&mut self,
vertex: Self::Vertex
) -> Result<FixedVertexHandle, InsertionError>
fn insert( &mut self, vertex: Self::Vertex ) -> Result<FixedVertexHandle, InsertionError>
source§fn fixed_undirected_edges(&self) -> FixedUndirectedEdgeIterator
fn fixed_undirected_edges(&self) -> FixedUndirectedEdgeIterator
source§fn fixed_directed_edges(&self) -> FixedDirectedEdgeIterator
fn fixed_directed_edges(&self) -> FixedDirectedEdgeIterator
source§fn fixed_all_faces(&self) -> FixedFaceIterator
fn fixed_all_faces(&self) -> FixedFaceIterator
source§fn fixed_inner_faces(&self) -> FixedInnerFaceIterator
fn fixed_inner_faces(&self) -> FixedInnerFaceIterator
source§fn all_vertices_on_line(&self) -> bool
fn all_vertices_on_line(&self) -> bool
true
if all vertices lie on a single line. Read more