Struct spade::ConstrainedDelaunayTriangulation
source · pub struct ConstrainedDelaunayTriangulation<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 constrained Delaunay triangulation.
A constrained Delaunay triangulation (CDT) is a triangulation that can contain constraint edges. These edges will always be present in the resulting triangulation.
Left: A CDT with 4 constraint edges. Right: The same triangulation without constraint edges
The resulting triangulation does not necessarily fulfill the Delaunay property.
This implementation currently supports only weakly intersecting constraints, thus, constraint edges are allowed to touch at their start or end point but are not allowed to intersect at any interior point.
The constrained triangulation shares most of the implementation of
the usual Delaunay triangulation, refer to DelaunayTriangulation
for more information about type parameters, iteration, performance
and more examples.
§Example
use spade::{ConstrainedDelaunayTriangulation, Point2, Triangulation};
let mut cdt = ConstrainedDelaunayTriangulation::<Point2<_>>::new();
let v0 = cdt.insert(Point2::new(0f64, 0.0))?;
let v1 = cdt.insert(Point2::new(1.0, 0.0))?;
cdt.add_constraint(v0, v1);
// Alternatively, consider using this shorthand
cdt.add_constraint_edge(Point2::new(1.0, 1.0), Point2::new(1.0, 0.0))?;
println!("Number of constraints: {}", cdt.num_constraints()); // 2 constraints
// Constraints are bidirectional!
assert!(cdt.exists_constraint(v1, v0));
assert!(cdt.exists_constraint(v0, v1));
// Check if a new constraint could be added
let from = Point2::new(1.0, -2.0);
let to = Point2::new(1.0, 0.0);
if !cdt.intersects_constraint(from, to) {
// No intersections, the edge can be added
cdt.add_constraint_edge(from, to)?;
}
§See also
Refer to Triangulation for most implemented methods on this type. Refer to DelaunayTriangulation for general information about using Delaunay triangulations.
Implementations§
source§impl<V, DE, UE, F, L> 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> ConstrainedDelaunayTriangulation<V, DE, UE, F, L>where
V: HasPosition,
DE: Default,
UE: Default,
F: Default,
L: HintGenerator<<V as HasPosition>::Scalar>,
sourcepub fn remove(&mut self, vertex: FixedVertexHandle) -> V
pub fn remove(&mut self, vertex: FixedVertexHandle) -> V
Removes a vertex from the triangulation.
This operation runs in O(n²), where n is the degree of the removed vertex.
§Handle invalidation
This method will invalidate all vertex, edge and face handles.
sourcepub fn num_constraints(&self) -> usize
pub fn num_constraints(&self) -> usize
Returns the number of constraint edges.
sourcepub fn is_constraint_edge(&self, edge: FixedUndirectedEdgeHandle) -> bool
pub fn is_constraint_edge(&self, edge: FixedUndirectedEdgeHandle) -> bool
Returns true
if a given edge is a constraint edge.
sourcepub fn exists_constraint(
&self,
from: FixedVertexHandle,
to: FixedVertexHandle
) -> bool
pub fn exists_constraint( &self, from: FixedVertexHandle, to: FixedVertexHandle ) -> bool
Checks if two vertices are connected by a constraint edge.
sourcepub fn can_add_constraint(
&self,
from: FixedVertexHandle,
to: FixedVertexHandle
) -> bool
pub fn can_add_constraint( &self, from: FixedVertexHandle, to: FixedVertexHandle ) -> bool
Checks if a constraint edge can be added.
Returns false
if the line from from
to to
intersects another
constraint edge.
sourcepub fn intersects_constraint(
&self,
line_from: Point2<V::Scalar>,
line_to: Point2<V::Scalar>
) -> bool
pub fn intersects_constraint( &self, line_from: Point2<V::Scalar>, line_to: Point2<V::Scalar> ) -> bool
Checks if a line intersects a constraint edge.
Returns true
if the edge from from
to to
intersects a
constraint edge.
sourcepub fn add_constraint_edges(
&mut self,
vertices: impl IntoIterator<Item = V>,
closed: bool
) -> Result<(), InsertionError>
pub fn add_constraint_edges( &mut self, vertices: impl IntoIterator<Item = V>, closed: bool ) -> Result<(), InsertionError>
Creates a several constraint edges by taking and connecting vertices from an iterator.
Every two sequential vertices in the input iterator will be connected by a constraint edge.
If closed
is set to true, the first and last vertex will also be connected.
§Special cases:
- Does nothing if input iterator is empty
- Only inserts the single vertex if the input iterator contains exactly one element
§Example
use spade::{ConstrainedDelaunayTriangulation, Point2};
const NUM_VERTICES: usize = 51;
let mut cdt = ConstrainedDelaunayTriangulation::<_>::default();
// Iterates through vertices on a circle
let vertices = (0..NUM_VERTICES).map(|i| {
let angle = std::f64::consts::PI * 2.0 * i as f64 / NUM_VERTICES as f64;
let (sin, cos) = angle.sin_cos();
Point2::new(sin, cos)
});
cdt.add_constraint_edges(vertices, true)?;
§Panics
Panics if any of the generated constraints intersects with any other constraint edge.
sourcepub fn add_constraint_edge(
&mut self,
from: V,
to: V
) -> Result<bool, InsertionError>
pub fn add_constraint_edge( &mut self, from: V, to: V ) -> Result<bool, InsertionError>
Insert two points and creates a constraint between them.
Returns true
if at least one constraint edge was added.
§Panics
Panics if the new constraint edge intersects with an existing constraint edge. Use can_add_constraint to check.
sourcepub fn add_constraint(
&mut self,
from: FixedVertexHandle,
to: FixedVertexHandle
) -> bool
pub fn add_constraint( &mut self, from: FixedVertexHandle, to: FixedVertexHandle ) -> bool
Adds a constraint edge between to vertices.
Returns true
if at least one constraint edge was added.
Note that the given constraint might be splitted into smaller edges
if a vertex in the triangulation lies exactly on the constraint edge.
Thus, cdt.exists_constraint(from, to)
is not necessarily true
after a call to this function.
Returns false and does nothing if from == to
.
§Panics
Panics if the new constraint edge intersects an existing constraint edge.
source§impl<V, DE, UE, F, L> ConstrainedDelaunayTriangulation<V, DE, UE, F, L>where
V: HasPosition + From<Point2<<V as HasPosition>::Scalar>>,
DE: Default,
UE: Default,
F: Default,
L: HintGenerator<<V as HasPosition>::Scalar>,
<V as HasPosition>::Scalar: Float,
impl<V, DE, UE, F, L> ConstrainedDelaunayTriangulation<V, DE, UE, F, L>where
V: HasPosition + From<Point2<<V as HasPosition>::Scalar>>,
DE: Default,
UE: Default,
F: Default,
L: HintGenerator<<V as HasPosition>::Scalar>,
<V as HasPosition>::Scalar: Float,
sourcepub fn refine(
&mut self,
parameters: RefinementParameters<V::Scalar>
) -> RefinementResult
pub fn refine( &mut self, parameters: RefinementParameters<V::Scalar> ) -> RefinementResult
Refines a triangulation by inserting additional points to improve the quality of its mesh.
Mesh quality, when applied to constrained delaunay triangulations (CDT), usually refers to how skewed its triangles are. A skewed triangle is a triangle with very large or very small (acute) inner angles. Some applications (e.g. interpolation and finite element methods) perform poorly in the presence of skewed triangles.
Refining by inserting additional points (called “steiner points”) may increase the minimum angle. The given RefinementParameters should be used to modify the refinement behavior.
§General usage
The vertex type must implement From<Point2<...>>
- otherwise, Spade cannot construct new steiner points at a
certain location. The refinement itself happens in place and will result in a valid CDT.
use spade::{ConstrainedDelaunayTriangulation, RefinementParameters, Point2, InsertionError, Triangulation};
fn get_refined_triangulation(vertices: Vec<Point2<f64>>) ->
Result<ConstrainedDelaunayTriangulation<Point2<f64>>, InsertionError>
{
let mut cdt = ConstrainedDelaunayTriangulation::bulk_load(vertices)?;
let result = cdt.refine(RefinementParameters::default());
if !result.refinement_complete {
panic!("Refinement failed - I should consider using different parameters.")
}
return Ok(cdt)
}
§Example image
Unrefined | Refined |
---|---|
A refinement example. The CDT on the left has some acute angles and skewed triangles. The refined CDT on the right contains several additional points that prevents such triangles from appearing while keeping all input vertices and constraint edges.
§Quality guarantees
Refinement will ensure that the resulting triangulation fulfills a few properties:
- The triangulation’s minimum angle will be larger than the angle specified by
with_angle_limit.
Exception: Acute input angles (small angles between initial constraint edges) cannot be maximized as the constraint edges must be kept intact. The algorithm will, for the most part, leave those places unchanged.
Exception: The refinement will often not be able to increase the minimal angle much above 30 degrees as every newly inserted steiner point may create additional skewed triangles. - The refinement will fullfil the Delaunay Property: Every triangle’s circumcenter will not contain any other vertex.
Exception: Refining with keep_constraint_edges cannot restore the Delaunay property if doing so would require splitting a constraint edge.
Exception: Refining with exclude_outer_faces will not restore the Delaunay property of any outer face. - Spade allows to specify a maximum allowed triangle area. The algorithm will attempt to subdivide any triangle with an area larger than this, independent of its smallest angle.
- Spade allows to specify a minimum required triangle area. The refinement will attempt to ignore any triangle with an area smaller than this parameter. This can prevent the refinement algorithm from over-refining in some cases.
§General limitations
The algorithm may fail to terminate in some cases for a minimum angle limit larger than 30 degrees. Such a limit can result in an endless loop: Every additionally inserted point creates more triangles that need to be refined.
To prevent this, spade limits the number of additionally inserted steiner points (see RefinementParameters::with_max_additional_vertices). However, this may leave the refinement in an incomplete state - some areas of the input mesh may not have been triangulated at all, some will be overly refined. Use RefinementResult::refinement_complete to identify if a refinement operation has succeeded without running out of vertices.
For mitigation, consider either lowering the minimum angle limit (see RefinementParameters::with_angle_limit) or introduce a minimum required area.
Meshes with very small input angles (angles between two constraint edges) may lead to poorly refined results. Please consider providing a bug report if you encounter an input mesh which you think isn’t refined well.
§Stability guarantees
While changing the interface of this method is considered to be a breaking change, changes to the specific refinement process (e.g. which faces are split in which order) are not. Any patch release may change how the same input mesh is being refined.
§References
This is an adaption of the classical refinement algorithms introduced by Jim Ruppert and Paul Chew.
For a good introduction to the topic, refer to the slides from a short course at the thirteenth and fourteenth International Meshing Roundtables (2005) by Jonathan Richard Shewchuk: https://people.eecs.berkeley.edu/~jrs/papers/imrtalk.pdf
Wikipedia: https://en.wikipedia.org/wiki/Delaunay_refinement
Trait Implementations§
source§impl<V, DE, UE, F, L> Clone for ConstrainedDelaunayTriangulation<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 ConstrainedDelaunayTriangulation<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) -> ConstrainedDelaunayTriangulation<V, DE, UE, F, L>
fn clone(&self) -> ConstrainedDelaunayTriangulation<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> Default 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> Default for ConstrainedDelaunayTriangulation<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 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> Triangulation for ConstrainedDelaunayTriangulation<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 = CdtEdge<UE>
type UndirectedEdge = CdtEdge<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