spade

Struct 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.

CDT No CDT

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>,

Source

pub fn bulk_load_cdt( vertices: Vec<V>, edges: Vec<[usize; 2]>, ) -> Result<Self, InsertionError>

Efficient bulk loading of a constraint delaunay triangulation, including both vertices and constraint edges.

The edges are given as pairs of vertex indices.

Note that the vertex order is not preserved by this function - iterating through all vertices will not result in the same sequence as the input vertices. Use ConstrainedDelaunayTriangulation::bulk_load_cdt_stable for a slower but order preserving variant.

Input vertices may have the same position. However, only one vertex for each position will be kept. Edges that go to a discarded vertex are rerouted and still inserted. It is arbitrary which duplicated vertex remains.

§Example
use spade::{ConstrainedDelaunayTriangulation, Point2, Triangulation};
let mut vertices = vec![
    Point2::new(0.0, 1.0),
    Point2::new(1.0, 2.0),
    Point2::new(3.0, -3.0),
    Point2::new(-1.0, -2.0),
    Point2::new(-4.0, -5.0),
];
let mut edges = vec![[0, 1], [1, 2], [2, 3], [3, 4]];
let cdt = ConstrainedDelaunayTriangulation::<_>::bulk_load_cdt(vertices.clone(), edges)?;

assert_eq!(cdt.num_vertices(), 5);
assert_eq!(cdt.num_constraints(), 4);
// The order will usually change
assert_ne!(cdt.vertices().map(|v| v.position()).collect::<Vec<_>>(), vertices);
§Panics

Panics if any constraint edges overlap. Panics if the edges contain an invalid index (out of range).

Source

pub fn bulk_load_cdt_stable( vertices: Vec<V>, edges: Vec<[usize; 2]>, ) -> Result<Self, InsertionError>

Stable bulk load variant that preserves the input vertex order

The resulting vertex set will be equal to the input vertex set if their positions are all distinct. See ConstrainedDelaunayTriangulation::bulk_load_cdt for additional details like panic behavior and duplicate handling.

§Example
use spade::{ConstrainedDelaunayTriangulation, Point2, Triangulation};
let mut vertices = vec![
    Point2::new(0.0, 1.0),
    Point2::new(1.0, 2.0),
    Point2::new(3.0, -3.0),
    Point2::new(-1.0, -2.0),
    Point2::new(-4.0, -5.0),
];
let mut edges = vec![[0, 1], [1, 2], [2, 3], [3, 4]];
let cdt = ConstrainedDelaunayTriangulation::<_>::bulk_load_cdt_stable(vertices.clone(), edges)?;

// The ordered will be preserved:
assert_eq!(cdt.vertices().map(|v| v.position()).collect::<Vec<_>>(), vertices);

It is fine to include vertex positions multiple times. The resulting order will be the same as if the duplicates were removed prior to insertion. However, it is unclear which duplicates are removed - e.g. do not assume that always the first duplicated vertex remains.

use spade::{ConstrainedDelaunayTriangulation, Point2, Triangulation};
let mut vertices = vec![
    Point2::new(0.0, 1.0),
    Point2::new(1.0, 2.0), // Duplicate
    Point2::new(1.0, 2.0),
    Point2::new(3.0, -3.0),
    Point2::new(3.0, -3.0), // Duplicate
    Point2::new(-4.0, -5.0),
];
let mut edges = vec![[0, 1], [2, 3], [4, 5]];
let cdt = ConstrainedDelaunayTriangulation::<_>::bulk_load_cdt_stable(vertices.clone(), edges)?;

// The choice of deduplicated vertices is arbitrary. In this example, dedup[1] and dedup[2] could
// have been swapped
let dedup = [
    Point2::new(0.0, 1.0),
    Point2::new(1.0, 2.0),
    Point2::new(3.0, -3.0),
    Point2::new(-4.0, -5.0),
];
assert_eq!(cdt.vertices().map(|v| v.position()).collect::<Vec<_>>(), dedup);
Source

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.

Source

pub fn num_constraints(&self) -> usize

Returns the number of constraint edges.

Source

pub fn is_constraint_edge(&self, edge: FixedUndirectedEdgeHandle) -> bool

Returns true if a given edge is a constraint edge.

Source

pub fn exists_constraint( &self, from: FixedVertexHandle, to: FixedVertexHandle, ) -> bool

Checks if two vertices are connected by a constraint edge.

Source

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.

Source

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.

Source

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.

Source

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.

Source

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 split 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. Use Self::try_add_constraint or Self::add_constraint_and_split to work around that.

Source

pub fn get_conflicting_edges_between_points( &self, from: Point2<<V as HasPosition>::Scalar>, to: Point2<<V as HasPosition>::Scalar>, ) -> impl Iterator<Item = DirectedEdgeHandle<'_, V, DE, CdtEdge<UE>, F>>

Returns all constraint edges that would prevent creating a new constraint between two points.

§See also

See also Self::get_conflicting_edges_between_vertices

Source

pub fn get_conflicting_edges_between_vertices( &self, from: FixedVertexHandle, to: FixedVertexHandle, ) -> impl Iterator<Item = DirectedEdgeHandle<'_, V, DE, CdtEdge<UE>, F>>

Returns all constraint edges that would prevent inserting a new constraint connecting two existing vertices.

§See also

See also Self::get_conflicting_edges_between_points

Source

pub fn remove_constraint_edge( &mut self, edge: FixedUndirectedEdgeHandle, ) -> bool

Removes a constraint edge.

Does nothing and returns false if the given edge is not a constraint edge. Otherwise, the edge is unmarked and the Delaunay property is restored in its vicinity.

Source

pub fn try_add_constraint( &mut self, from: FixedVertexHandle, to: FixedVertexHandle, ) -> Vec<FixedDirectedEdgeHandle>

Attempts to add a constraint edge. Leaves the triangulation unchanged if the new edge would intersect an existing constraint edge.

Returns all constraint edges that connect from and to. This includes any constraint edge that was already present. Multiple edges are returned if the line from from to to intersects an existing vertex. Returns an empty list if the new constraint would intersect any existing constraint or if from == to.

§Example
use spade::{ConstrainedDelaunayTriangulation, Point2, Triangulation};
let mut cdt = ConstrainedDelaunayTriangulation::<Point2<_>>::new();
let v0 = cdt.insert(Point2::new(-1.0, 0.0))?;
let v1 = cdt.insert(Point2::new(1.0, 0.0))?;
let v2 = cdt.insert(Point2::new(0.0, 1.0))?;
let v3 = cdt.insert(Point2::new(0.0, -1.0))?;
let first_constraints = cdt.try_add_constraint(v2, v3);
let second_constraints = cdt.try_add_constraint(v0, v1);

// The first constraint edge can be added as there are no intersecting constraint edges
assert_eq!(first_constraints.len(), 1);
let edge = cdt.directed_edge(first_constraints[0]);
assert_eq!(edge.from().fix(), v2);
assert_eq!(edge.to().fix(), v3);

// The second edge should not be created as it intersects the first edge.
assert!(second_constraints.is_empty());

// Consider comparing this to the number of constraints prior to calling
// `try_add_constraint` to check if any new constraint edge was created.
assert_eq!(cdt.num_constraints(), 1);
Source§

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

Source

pub fn add_constraint_and_split<C>( &mut self, from: FixedVertexHandle, to: FixedVertexHandle, vertex_constructor: C, ) -> Vec<FixedDirectedEdgeHandle>
where C: Fn(Point2<<V as HasPosition>::Scalar>) -> V,

Adds a constraint to the triangulation. Splits any existing constraint edge that would intersect the new constraint edge.

The vertex_constructor closure is used to convert the position of the intersection into a vertex. The returned vertex must have exactly the same position as the argument of the closure.

Returns all constraint edges that connect from and to. This includes any constraint edge that was already present. Multiple edges are returned if the line from from to to intersects any existing vertex or any existing constraint edge. Returns an empty list if from == to.

§Image example

This is an input CDT with 3 constraints:

0 1

Calling add_constraint_and_split(v0, v1, ...) will result in this CDT:

0 1
§Code example
 use spade::{ConstrainedDelaunayTriangulation, Point2, Triangulation};
 use spade::handles::FixedVertexHandle;
 let mut cdt = ConstrainedDelaunayTriangulation::<Point2<_>>::new();
 let v0 = cdt.insert(Point2::new(-1.0, 0.0))?;
 let v1 = cdt.insert(Point2::new(1.0, 0.0))?;
 let v2 = cdt.insert(Point2::new(0.0, 1.0))?;
 let v3 = cdt.insert(Point2::new(0.0, -1.0))?;
 cdt.add_constraint(v2, v3);

 // Should create a new split vertex at the origin
 let second_constraints = cdt.add_constraint_and_split(v0, v1, |v| v);

 // Expect one additional point introduced by splitting the first constraint edge:
 assert_eq!(cdt.num_vertices(), 5);

 let v4 = FixedVertexHandle::from_index(4); // Newly created

 // Expect 4 constraints as the first constraint was split:
 assert_eq!(cdt.num_constraints(), 4);

 // The second edge should consist of two edges, v0 -> v4 and v4 -> v1
 assert_eq!(second_constraints.len(), 2);

 let [e0, e1] = [second_constraints[0], second_constraints[1]];
 let [e0, e1] = [e0, e1].map(|e| cdt.directed_edge(e));

 assert_eq!(e0.from().fix(), v0);
 assert_eq!(e0.to().fix(), v4);
 assert_eq!(e1.from().fix(), v4);
 assert_eq!(e1.to().fix(), v1);
§Precision warning

Intersection points may not exactly lie on the line between from and to, either due to loss of precision or as the exact value may not be representable with the underlying floating point number.

Thus, iterating a LineIntersectionIterator::new_from_handles(&cdt, from, to) will often not return only Intersection::EdgeOverlap as would be expected. Instead, use the returned Vec to identify the edges that form the new constraint. The absolute deviation from the correct position should be small, especially when using f64 coordinates as storage type.

Source§

impl<V, DE, UE, F, L> ConstrainedDelaunayTriangulation<V, DE, UE, F, L>

Source

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
UnrefinedRefined

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,

Source§

fn clone(&self) -> ConstrainedDelaunayTriangulation<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> Default for ConstrainedDelaunayTriangulation<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 ConstrainedDelaunayTriangulation<V, DE, UE, F, L>
where V: HasPosition, DE: Default, UE: Default, F: Default, L: HintGenerator<<V as HasPosition>::Scalar>,

Source§

type Vertex = V

The triangulation’s vertex type.
Source§

type DirectedEdge = DE

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

type UndirectedEdge = CdtEdge<UE>

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

type Face = F

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

type HintGenerator = L

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

fn clear(&mut self)

Removes all edges, faces and vertices from the triangulation. Read more
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 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 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 get_vertex( &self, handle: FixedVertexHandle, ) -> Option<VertexHandle<'_, Self::Vertex, Self::DirectedEdge, Self::UndirectedEdge, Self::Face>>

Get a vertex by its index 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. Read more
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> Freeze for ConstrainedDelaunayTriangulation<V, DE, UE, F, L>
where L: Freeze,

§

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

§

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

§

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

§

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

§

impl<V, DE, UE, F, L> UnwindSafe for ConstrainedDelaunayTriangulation<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> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dst: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dst. 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,

Source§

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>,

Source§

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>,

Source§

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.