Algorithm Internals

Here we give some functions that are used for some of the algorithms in this package.

Unconstrained Triangulations

Here are some of the internal functions used for the Bowyer-Watson algorithm, i.e. for the computation of an unconstrained triangulation.

DelaunayTriangulation.postprocess_triangulate!Function
postprocess_triangulate!(tri; delete_ghosts=false, delete_empty_features=true, recompute_representative_points=true)

Postprocesses the triangulation tri after it has been constructed using triangulate. This includes:

source
DelaunayTriangulation.add_point_bowyer_watson!Method
add_point_bowyer_watson!(tri::Triangulation, new_point, initial_search_point::I, rng::Random.AbstractRNG=Random.default_rng(), update_representative_point=true, store_event_history=Val(false), event_history=nothing, peek::P=Val(false), predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Triangle

Adds new_point into tri.

Arguments

  • tri: The triangulation.
  • new_point: The point to insert.
  • initial_search_point::I: The vertex to start the point location with find_triangle at. See get_initial_search_point.
  • rng::Random.AbstractRNG: The random number generator to use.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • update_representative_point=true: If true, then the representative point is updated. See update_centroid_after_addition!.
  • store_event_history=Val(false): If true, then the event history from the insertion is stored.
  • event_history=nothing: The event history to store the event history in. Should be an InsertionEventHistory if store_event_history is true, and false otherwise.
  • peek=Val(false): Whether to actually add new_point into tri, or just record into event_history all the changes that would occur from its insertion.

Output

  • V: The triangle in tri containing new_point.

Extended help

This function works as follows:

  1. First, the triangle containing the new point, V, is found using find_triangle.
  2. Once the triangle is found, we call into add_point_bowyer_watson_and_process_after_found_triangle to properly insert the point.
  3. Inside add_point_bowyer_watson_and_process_after_found_triangle, we first call into add_point_bowyer_watson_after_found_triangle to add the point into the cavity. We then call into add_point_bowyer_watson_onto_segment to make any changes necessary incase the triangulation is constrained and new_point lies on a segment, since the depth-first search of the triangles containing new_point in its circumcenter must be performed on each side of the segment that new_point lies on.

The function add_point_bowyer_watson_dig_cavities! is the main workhorse of this function from add_point_bowyer_watson_after_found_triangle. See its docstring for the details.

source
DelaunayTriangulation.add_point_bowyer_watson_dig_cavities!Method
add_point_bowyer_watson_dig_cavities!(tri::Triangulation, new_point::N, V, q, flag, update_representative_point=true, store_event_history=Val(false), event_history=nothing, peek::F=Val(false), predicates::AbstractPredicateKernel=AdaptiveKernel()) where {N,F}

Deletes all the triangles in tri whose circumcircle contains new_point. This leaves behind a polygonal cavity, whose boundary edges are then connected to new_point, restoring the Delaunay property from new_point's insertion.

Arguments

  • tri: The Triangulation.
  • new_point::N: The point to insert.
  • V: The triangle in tri containing new_point.
  • q: The point to insert.
  • flag: The position of q relative to V. See point_position_relative_to_triangle.
  • update_representative_point=true: If true, then the representative point is updated. See update_centroid_after_addition!.
  • store_event_history=Val(false): If true, then the event history from the insertion is stored.
  • event_history=nothing: The event history to store the event history in. Should be an InsertionEventHistory if store_event_history is true, and false otherwise.
  • peek=Val(false): Whether to actually add new_point into tri, or just record into event_history all the changes that would occur from its insertion.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

There are no changes, but tri is updated in-place.

Extended help

This function works as follows:

  1. To dig the cavity, we call dig_cavity! on each edge of V, stepping towards the adjacent triangles to excavate the cavity recursively.
  2. Once the cavity has been excavated, extra care is needed in case is_on(flag), meaning new_point is on one of the edges of V. In particular, extra care is needed if: is_on(flag) && (is_boundary_triangle(tri, V) || is_ghost_triangle(V) && !is_boundary_node(tri, new_point)[1]). The need for this check is in case new_point is on a boundary edge already exists, since we need to fix the associated ghost edges. For example, a boundary edge (i, k) might have been split into the edges (i, j) and (j, k), which requires that the ghost triangle (k, i, g) be split into (j, i, g) and (k, j, g), where g is the ghost vertex. This part of the function will fix this case. The need for is_ghost_triangle(V) && !is_boundary_node(tri, new_point)[1] is in case the ghost edges were already correctly added. Nothing happens if the edge of V that new_point is on is not the boundary edge.
source
DelaunayTriangulation.dig_cavity!Method
dig_cavity!(tri::Triangulation, r, i, j, ℓ, flag, V, store_event_history=Val(false), event_history=nothing, peek::F=Val(false), predicates::AbstractPredicateKernel=AdaptiveKernel()) where {F}

Excavates the cavity in tri through the edge (i, j), stepping towards the adjacent triangles to excavate the cavity recursively, eliminating all triangles containing r in their circumcircle.

Arguments

  • tri: The Triangulation.
  • r: The new point being inserted.
  • i: The first vertex of the edge (i, j).
  • j: The second vertex of the edge (i, j).
  • : The vertex adjacent to (j, i), so that the triangle being stepped into is (j, i, ℓ).
  • flag: The position of r relative to V. See point_position_relative_to_triangle.
  • V: The triangle in tri containing r.
  • store_event_history=Val(false): If true, then the event history from the insertion is stored.
  • event_history=nothing: The event history to store the event history in. Should be an InsertionEventHistory if store_event_history is true, and false otherwise.
  • peek=Val(false): Whether to actually add new_point into tri, or just record into event_history all the changes that would occur from its insertion.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

There are no changes, but tri is updated in-place.

Extended help

This function works as follows:

  1. First, we check if is 0, which would imply that the triangle (j, i, ℓ) doesn't exist as it has already been deleted. If this is the check, we return and stop digging.
  2. If the edge (i, j) is not a segment, r is inside of the circumcenter of (j, i, ℓ), and is not a ghost vertex, then we can step forward into the next triangles. In particular, we delete the triangle (j, i, ℓ) from tri and then call dig_cavity! again on two edges of (j, i, ℓ) other than (j, i).
  3. If we did not do step 2, then this means that we are on the edge of the polygonal cavity; this also covers the case that is a ghost vertex, i.e. we are at the boundary of the triangulation. There are two cases to consider here. Firstly, if (i, j) is not a segment, we just need to add the triangle (r, i, j) into tri to connect r to this edge of the polygonal cavity. If (i, j) is a segment, then the situation is more complicated. In particular, r being on an edge of V might imply that we are going to add a degenerate triangle (r, i, j) into tri, and so this needs to be avoided. So, we check if is_on(flag) && contains_segment(tri, i, j) and, if the edge that r is on is (i, j), we add the triangle (r, i, j). Otherwise, we do nothing.
source
DelaunayTriangulation.enter_cavityFunction
enter_cavity(tri::Triangulation, r, i, j, ℓ, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Bool

Determines whether to enter the cavity in tri through the edge (i, j) when inserting r into the triangulation.

Arguments

  • tri: The Triangulation.
  • r: The new point being inserted.
  • i: The first vertex of the edge (i, j).
  • j: The second vertex of the edge (i, j).
  • : The vertex adjacent to (j, i), so that the triangle being stepped into is (j, i, ℓ).
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

source
DelaunayTriangulation.get_initial_search_pointMethod
get_initial_search_point(tri::Triangulation, num_points, new_point, insertion_order, num_sample_rule::F, rng, try_last_inserted_point) where {F} -> Vertex

For a given iteration of the Bowyer-Watson algorithm, finds the point to start the point location with find_triangle at.

Arguments

  • tri: The triangulation.
  • num_points: The number of points currently in the triangulation.
  • new_point: The point to insert.
  • insertion_order: The insertion order of the points. See get_insertion_order.
  • num_sample_rule::F: The rule to use to determine the number of points to sample. See default_num_samples for the default.
  • rng::Random.AbstractRNG: The random number generator to use.
  • try_last_inserted_point: If true, then the last inserted point is also considered as the start point.

Output

  • initial_search_point: The vertex to start the point location with find_triangle at.
source
DelaunayTriangulation.get_initial_triangleFunction
get_initial_triangle(tri::Triangulation, insertion_order, itr=0) -> Triangle

Gets the initial triangle for the Bowyer-Watson algorithm.

Arguments

  • tri: The triangulation.
  • insertion_order: The insertion order of the points. See get_insertion_order.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • itr=0: To avoid issues with degenerate triangles and infinite loops, this counts the number of times insertion_order had to be shifted using circshift! to find an initial non-degenerate triangle.

Output

  • initial_triangle: The initial triangle.
source
DelaunayTriangulation.get_insertion_orderMethod
get_insertion_order(points, randomise, skip_points, ::Type{I}, rng) where {I} -> Vector{I}
get_insertion_order(tri::Triangulation, randomise, skip_points, rng) -> Vector{I}

Gets the insertion order for points into a triangulation.

Arguments

  • points: The points to insert.
  • randomise: If true, then the insertion order is randomised. Otherwise, the insertion order is the same as the order of the points.
  • skip_points: The points to skip.
  • I::Type{I}: The type of the vertices.
  • rng::Random.AbstractRNG: The random number generator to use.

Output

  • order: The order to insert the points in.
Mutation of `order`

This order might be mutated (by circshift!) in get_initial_triangle.

source
DelaunayTriangulation.initialise_bowyer_watson!Function
initialise_bowyer_watson!(tri::Triangulation, insertion_order, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Triangulation

Initialises the Bowyer-Watson algorithm.

Arguments

  • tri: The triangulation.
  • insertion_order: The insertion order of the points. See get_insertion_order.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

tri is updated in place to contain the initial triangle from which the Bowyer-Watson algorithm starts.

source
DelaunayTriangulation.unconstrained_triangulation!Method
unconstrained_triangulation!(tri::Triangulation; kwargs...)

Computes the unconstrained Delaunay triangulation of the points in tri.

Arguments

  • tri: The triangulation.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • randomise=true: If true, then the insertion order is randomised. Otherwise, the insertion order is the same as the order of the points.
  • skip_points=(): The vertices to skip.
  • num_sample_rule::M=default_num_samples: The rule to use to determine the number of points to sample. See default_num_samples for the default.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use.
  • insertion_order=get_insertion_order(tri, randomise, skip_points, rng): The insertion order of the points. See get_insertion_order.

Outputs

There is no output, but tri is updated in-place.

source

Triangulation Rectangular Domains

Here are some of the internal functions used for triangulate_rectangle.

DelaunayTriangulation.get_lattice_trianglesFunction
get_lattice_triangles(nx, ny, Ts, V)

Computes the triangles defining a lattice with nx and ny points in the x- and y-directions, respectively.

See triangulate_rectangle.

Arguments

  • nx: The number of x points in the lattice.
  • ny: The number of y points in the lattice.
  • Ts: The type to use for representing a collection of triangles.
  • V: The type to use for representing an individual triangle.

Outputs

  • T: The collection of triangles.
  • sub2ind: A map that takes cartesian indices (i, j) into the associated linear index along the lattice. See LinearIndices.
source
DelaunayTriangulation.get_lattice_pointsFunction
get_lattice_points(a, b, c, d, nx, ny, sub2ind)

Returns the points on a lattice [a, b] × [c, d].

See triangulate_rectangle.

Arguments

  • a: The minimum x-coordinate.
  • b: The maximum x-coordinate.
  • c: The minimum y-coordinate.
  • d: The maximum y-coordinate.
  • nx: The number of points in the x-direction.
  • ny: The number of points in the y-direction.
  • sub2ind: The map returned from get_lattice_triangles.

Outputs

  • points: The points on the lattice, where points[sub2ind[CartesianIndex(i, j)]] is the point at the ith x point and the jth y point,
source
DelaunayTriangulation.get_lattice_boundaryFunction
get_lattice_boundary(nx, ny, sub2ind, single_boundary, IntegerType)

Returns the boundary nodes defining a lattice.

See triangulate_rectangle.

Arguments

  • nx: The number of x points in the lattice.
  • ny: The number of y points in the lattice.
  • sub2ind: The map returned from get_lattice_triangles.
  • single_boundary: If true, then the boundary nodes are stored as a contiguous section. Otherwise, the boundary is split into four sections, in the order bottom, right, top, left.
  • IntegerType: The type used for representing vertices.

Outputs

  • boundary_nodes: The boundary nodes, returned according to single_boundary as described above.
source

Triangulating Convex Polygons

Here are some functions used in triangulate_convex for triangulating a convex polygon.

DelaunayTriangulation.add_point_convex_triangulation!Function
add_point_convex_triangulation!(tri::Triangulation, u, v, w, S, predicates::AbstractPredicateKernel=AdaptiveKernel())

Adds the point u into the triangulation tri.

Arguments

  • tri::Triangulation: The Triangulation.
  • u: The vertex to add.
  • v: The vertex next to u.
  • w: The vertex previous to u.
  • S: The set of vertices of the polygon.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, as tri is modified in-place.

Extended help

This function forms part of Chew's algorithm for triangulating a convex polygon. There are some important points to make.

  1. Firstly, checking that x = get_adjacent(tri, w, v) is needed to prevent the algorithm from exiting the polygon.

This is important in case this algorithm is used as part of delete_point!. When you are just triangulating a convex polygon by itself, this checked is the same as checking edge_exists(tri, w, v).

  1. For this method to be efficient, the set x ∈ S must be O(1) time. This is why we use a Set type for S.
  2. The algorithm is recursive, recursively digging further through the polygon to find non-Delaunay edges to adjoins with u.
source
DelaunayTriangulation.delete_vertices_in_random_order!Method
delete_vertices_in_random_order!(list::Triangulation, tri::ShuffledPolygonLinkedList, rng, predicates::AbstractPredicateKernel)

Deletes the vertices of the polygon represented by list in random order, done by switching the pointers of the linked list. Only three vertices will survive. If these these three vertices are collinear, then the deletion is attempted again after reshuffling the vertices.

Arguments

  • list::ShuffledPolygonLinkedList: The linked list representing the polygon to be deleted.
  • tri::Triangulation: The Triangulation.
  • rng::Random.AbstractRNG: The random number generator used to shuffle the vertices of S.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, as list is modified in-place.

source
DelaunayTriangulation.postprocess_triangulate_convex!Method
postprocess_triangulate_convex!(tri::Triangulation, S; delete_ghosts, delete_empty_features)

Postprocesses the completed triangulation tri of the convex polygon S.

Arguments

Keyword Arguments

  • delete_ghosts=false: If true, the ghost triangles are deleted after triangulation.
  • delete_empty_features=true: If true, the empty features are deleted after triangulation.

Output

There are no output, as tri is modified in-place. The postprocessing that is done is:

  1. The convex hull of tri is set to S.
  2. The ghost triangles are deleted if delete_ghosts=true.
  3. The empty features are deleted if delete_empty_features=true.
  4. The representative points are set to the centroid of S.
source
DelaunayTriangulation.triangulate_convex!Method
triangulate_convex!(tri::Triangulation, S; rng::Random.AbstractRNG=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel())

Triangulates the convex polygon S in-place into tri.

Arguments

  • tri::Triangulation: The triangulation to be modified.
  • S: A convex polygon represented as a vector of vertices. The vertices should be given in counter-clockwise order, and must not be circular so that S[begin] ≠ S[end].

Keyword Arguments

  • rng=Random.default_rng(): The random number generator used to shuffle the vertices of S before triangulation.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, as tri is updated in-place. This function does not do any post-processing, e.g. deleting any ghost triangles. This is done by triangulate_convex or postprocess_triangulate_convex!.

source
DelaunayTriangulation.triangulate_convexMethod
triangulate_convex(points, S; delete_ghosts=false, delete_empty_features=true, rng=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel(), kwargs...) -> Triangulation

Triangulates the convex polygon S.

Arguments

  • points: The point set corresponding to the vertices in S.
  • S: A convex polygon represented as a vector of vertices. The vertices should be given in counter-clockwise order, and must not be circular so that S[begin] ≠ S[end].

Keyword Arguments

  • delete_ghosts=false: If true, the ghost triangles are deleted after triangulation.
  • delete_empty_features=true: If true, the empty features are deleted after triangulation.
  • rng=Random.default_rng(): The random number generator used to shuffle the vertices of S before triangulation.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • kwargs...: Additional keyword arguments passed to Triangulation.

Output

  • tri::Triangulation: The triangulated polygon.
source

Constrained Triangulations

There are many functions to list for the computation of a constrained triangulation.

DelaunayTriangulation.add_segment!Method
add_segment!(tri::Triangulation, segment; predicates::AbstractPredicateKernel=AdaptiveKernel(), rng::Random.AbstractRNG=Random.default_rng())
add_segment!(tri::Triangulation, i, j; predicates::AbstractPredicateKernel=AdaptiveKernel(), rng::Random.AbstractRNG=Random.default_rng())

Adds segment = (i, j) to tri.

Arguments

  • tri::Triangulation: The Triangulation.
  • segment: The segment to add. The second method uses (i, j) to represent the segment instead.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • rng::AbstractRNG=Random.default_rng(): The RNG object.

Outputs

There is no output, but tri will be updated so that it now contains segment.

source
DelaunayTriangulation.add_segment_to_list!Method
add_segment_to_list!(tri::Triangulation, e)

Adds e to get_interior_segments(tri) and get_all_segments(tri) if it, or reverse_edge(e), is not already in the sets.

Arguments

Outputs

There is no output, but tri will be updated so that e is in get_interior_segments(tri) and get_all_segments(tri).

source
DelaunayTriangulation.optimise_edge_orderMethod
optimise_edge_order(tri::Triangulation, segment) -> Edge

Optimises the orientation of segment for inserting it into the triangulation.

Arguments

  • tri::Triangulation: The Triangulation.
  • segment: The segment to arrange.

Outputs

  • e: If segment is a boundary edge, then e = segment, Otherwise, e = sort_edge_by_degree(tri, segment) so that initial(e) has the smaller degree of the two vertices.
source
DelaunayTriangulation.add_point_cavity_cdt!Function
add_point_cavity_cdt!(tri::Triangulation, u, v, w, marked_vertices)

Adds a point to the cavity V left behind when deleting triangles intersected in a triangulation by an edge, updating tri to do so.

Arguments

  • tri::Triangulation: The Triangulation to update.
  • u: The vertex to add.
  • v: The vertex along the polygon that is next to u.
  • w: The vertex along the polygon that is previous to u.
  • marked_vertices: Cache for marking vertices to re-triangulate during the triangulation. This gets mutated.

Outputs

There is no output, but tri is updated in-place, as is marked_vertices if necessary.

source
DelaunayTriangulation.connect_segments!Method
connect_segments!(segments)

Connects the ordered vector of segments so that the endpoints all connect, preserving order.

Example

julia> using DelaunayTriangulation

julia> segments = [(7, 12), (12, 17), (17, 22), (32, 37), (37, 42), (42, 47)];

julia> DelaunayTriangulation.connect_segments!(segments)
7-element Vector{Tuple{Int64, Int64}}:
 (7, 12)
 (12, 17)
 (17, 22)
 (22, 32)
 (32, 37)
 (37, 42)
 (42, 47)
source
DelaunayTriangulation.constrained_triangulation!Method
constrained_triangulation!(tri::Triangulation, segments, boundary_nodes, predicates::AbstractPredicateKernel, full_polygon_hierarchy; rng=Random.default_rng(), delete_holes=true) -> Triangulation

Creates a constrained triangulation from tri by adding segments and boundary_nodes to it. This will be in-place, but a new triangulation is returned to accommodate the changed types.

Arguments

  • tri::Triangulation: The Triangulation.
  • segments: The interior segments to add to the triangulation.
  • boundary_nodes: The boundary nodes to add to the triangulation.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • full_polygon_hierarchy: The PolygonHierarchy defining the boundary. This will get copied into the existing polygon hierarchy.

Keyword Arguments

  • rng=Random.default_rng(): The random number generator to use.
  • delete_holes=true: Whether to delete holes in the triangulation. See delete_holes!.

Outputs

source
DelaunayTriangulation.delete_polygon_vertices_in_random_order!Method
delete_polygon_vertices_in_random_order!(list::ShuffledPolygonLinkedList, tri::Triangulation, u, v, rng::Random.AbstractRNG, predicates::AbstractPredicateKernel)

Deletes vertices from the polygon defined by list in a random order.

Arguments

  • list::ShuffledPolygonLinkedList: The linked list of polygon vertices.
  • tri::Triangulation: The Triangulation.
  • u, v: The vertices of the segment (u, v) that was inserted in order to define the polygon V = list.S.
  • rng::Random.AbstractRNG: The random number generator to use.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, but list is updated in-place.

source
DelaunayTriangulation.extend_segments!Method
extend_segments!(segments, segment)

Given an ordered vector of segments, ensures that they also represent the replacement of segment. In particular, suppose segments represents the sequence of edges

    ---(i₁, i₂)---(i₂, i₃)---(⋯, ⋯)---(iₙ₋₁, iₙ)---

and segment is (i₀, iₙ₊₁). Then the extended sequence becomes

    ---(i₀, i₁)---(i₁, i₂)---(i₂, i₃)---(⋯, ⋯)---(iₙ₋₁, iₙ)---(iₙ, iₙ₊₁)---

Example

julia> using DelaunayTriangulation

julia> segments = [(2, 7), (7, 12), (12, 49)];

julia> segment = (1, 68);

julia> DelaunayTriangulation.extend_segments!(segments, segment)
5-element Vector{Tuple{Int64, Int64}}:
 (1, 2)
 (2, 7)
 (7, 12)
 (12, 49)
 (49, 68)
source
DelaunayTriangulation.fix_segments!Method
fix_segments!(segments, bad_indices)

Fixes the overlapping segments in segments, referred to via bad_indices, by connecting consecutive edges where needed.

Arguments

  • segments: The segments to fix.
  • bad_indices: The indices of the segments to fix.

Outputs

There are no outputs as segments is updated in-place.

Example

julia> using DelaunayTriangulation

julia> segments = [(2, 15), (2, 28), (2, 41)]; # the edges all start with 2, so they are not actual segments in the triangulation, and so must be fixed

julia> bad_indices = [1, 2, 3];

julia> DelaunayTriangulation.fix_segments!(segments, bad_indices)
3-element Vector{Tuple{Int64, Int64}}:
 (2, 15)
 (15, 28)
 (28, 41)

julia> segments = [(2, 7), (2, 12), (12, 17), (2, 22), (2, 27), (2, 32), (32, 37), (2, 42), (42, 47)];

julia> bad_indices = [2, 4, 5, 6, 8]
5-element Vector{Int64}:
 2
 4
 5
 6
 8

julia> DelaunayTriangulation.fix_segments!(segments, bad_indices)
9-element Vector{Tuple{Int64, Int64}}:
 (2, 7)
 (7, 12)
 (12, 17)
 (17, 22)
 (22, 27)
 (27, 32)
 (32, 37)
 (37, 42)
 (42, 47)
source
DelaunayTriangulation.is_vertex_closer_than_neighboursMethod
is_vertex_closer_than_neighbours([predicates::AbstractPredicateKernel=AdaptiveKernel(),] tri::Triangulation, u, v, jᵢ, jᵢ₋₁, jᵢ₊₁) -> Bool
is_vertex_closer_than_neighbours([predicates::AbstractPredicateKernel=AdaptiveKernel(),] tri::Triangulation, list::ShuffledPolygonLinkedList, u, v, j) -> Bool

Tests if the vertex jᵢ is closer to the line (u, v) than its neighbours jᵢ₋₁ and jᵢ₊₁, assuming all these vertices are to the left of the line.

See also point_closest_to_line.

Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • tri::Triangulation: The Triangulation.
  • u, v: The vertices of the line.
  • jᵢ, jᵢ₋₁, jᵢ₊₁: The vertices to compare.

The second method extracts these latter two vertices using the doubly-linked list of vertices.

Outputs

  • flag: Whether jᵢ is closer to the line than jᵢ₋₁ and jᵢ₊₁.
source
DelaunayTriangulation.locate_intersecting_trianglesFunction
locate_intersecting_triangles(tri::Triangulation, e, rotate=Val(true), rng::Random.AbstractRNG=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Vector, Vector, Vector, Vector)

Find all the triangles intersected by an edge e.

See also find_triangle.

Arguments

  • tri::Triangulation: The Triangulation.
  • e: The edge going through the triangulation.
  • rotate=Val(true): Whether to rotate the edge so that the minimum degree vertex of e is first.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • intersecting_triangles: The intersected triangles.
  • collinear_segments: Segments that are collinear with e.
  • left_vertices: The vertices of the intersected triangles that are left of e.
  • right_vertices: The vertices of the intersected triangles that are right of e.
source
DelaunayTriangulation.merge_segmentsMethod
merge_segments(tri::Triangulation, ghost_vertex_map) -> Edges

Creates a set of edges that merge all the boundary nodes in tri as well as its interior segments into a single collection.

Arguments

  • tri::Triangulation: The Triangulation.
  • ghost_vertex_map: The ghost vertex map to use.

Outputs

  • all_segments: The set of edges that merge all the boundary nodes in tri as well as its interior segments into a single collection, with type equal to that of get_interior_segments(tri)'s.
source
DelaunayTriangulation.prepare_vertex_linked_listMethod
prepare_vertex_linked_list(V) -> ShuffledPolygonLinkedList

Given a list of polygon vertices V, returns the doubly-linked list of polygon vertices.

Arguments

  • V: The list of polygon vertices.

Outputs

  • list::ShuffledPolygonLinkedList: A ShuffledPolygonLinkedList. In list, prev[begin], prev[end], next[begin], and next[end] are all 0 as are shuffled_indices[begin] and shuffled_indices[end]. Moreover, shuffled_indices will not have been shuffled yet.
source
DelaunayTriangulation.process_collinear_segments!Method
process_intersecting_triangles!(tri::Triangulation, e, collinear_segments; predicates::AbstractPredicateKernel=AdaptiveKernel(), rng::Random.AbstractRNG=Random.default_rng()) -> Bool

Given segments in collinear_segments that are collinear with an edge e, updates tri so that this edge e is instead split so that it is instead represented by collinear_segments. These new segments will be placed into the triangulation using add_segment!.

The predicates::AbstractPredicateKernel argument defines the method for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

See also connect_segments!, extend_segments!, split_segment! and split_boundary_edge_at_collinear_segments!.

source
DelaunayTriangulation.remake_triangulation_with_constraintsMethod
remake_triangulation_with_constraints(tri::Triangulation, segments, boundary_nodes) -> (Dict, Dict, Triangulation)

Remakes the triangulation tri so that it contains segments and boundary_nodes in its fields.

See also replace_ghost_vertex_information.

Arguments

  • tri::Triangulation: The triangulation to remake.
  • segments: The segments to add to the triangulation.
  • boundary_nodes: The boundary nodes to add to the triangulation.

Outputs

  • new_ghost_vertex_map: The new ghost vertex map. This will not yet be added to the triangulation.
  • new_ghost_vertex_ranges: The new ghost vertex ranges. This will not yet be added to the triangulation.
  • new_tri::Triangulation: The new triangulation, now containing boundary_nodes in the boundary_nodes field and segments in the interior_segments field.
source
DelaunayTriangulation.replace_ghost_vertex_informationMethod
replace_ghost_vertex_information(tri::Triangulation, ghost_vertex_map, ghost_vertex_ranges) -> Triangulation

Replaces the ghost vertex information in tri with ghost_vertex_map and ghost_vertex_ranges, using the results from remake_triangulation_with_constraints.

Arguments

  • tri::Triangulation: The triangulation to remake.
  • ghost_vertex_map: The ghost vertex map to add to the triangulation.
  • ghost_vertex_ranges: The ghost vertex ranges to add to the triangulation.

Outputs

  • new_tri::Triangulation: The new triangulation, now containing ghost_vertex_map in the ghost_vertex_map field and ghost_vertex_ranges in the ghost_vertex_ranges field.
source
DelaunayTriangulation.retriangulate_fan!Method
retriangulate_fan!(tri::Triangulation, tri_fan::Triangulation, fan, fan_triangles; predicates::AbstractPredicateKernel=AdaptiveKernel(), rng::Random.AbstractRNG=Random.default_rng())

Given a sorted set of vertices fan in a fan of triangles associated with fan_triangles, retriangulates the fan, updating tri to do so and using tri_fan as a temporary triangulation. (This implements Lines 17–19 and Line 28 of the algorithms in this paper.)

The predicates argument defines the method for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

source
DelaunayTriangulation.select_random_vertexMethod
select_random_vertex(tri::Triangulation, list::ShuffledPolygonLinkedList, u, v, range, rng) -> Vertex

Selects a random vertex that is not closer to the line (u, v) than both of its neighbours.

Arguments

  • tri::Triangulation: The Triangulation.
  • list::ShuffledPolygonLinkedList: The linked list of polygon vertices.
  • u, v: The vertices of the line.
  • range: The range of indices of the vertices to select from.
  • rng::Random.AbstractRNG: The random number generator to use.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • j: The selected vertex.
source
DelaunayTriangulation.setup_cavity_cdtMethod
setup_cavity_cdt(tri::Triangulation, V; rng::Random.AbstractRNG=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel()) -> ShuffledPolygonLinkedList

Prepares the linked list required for triangulating a cavity excavated by segment insertion in a constrained triangulation.

See also prepare_vertex_linked_list and delete_polygon_vertices_in_random_order!.

Arguments

  • tri::Triangulation: The Triangulation.
  • V: The list of polygon vertices, given as a counter-clockwise list of vertices, defining the cavity.

Keyword Arguments

  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • list::ShuffledPolygonLinkedList: The linked list of polygon vertices representing the cavity.
source
DelaunayTriangulation.sort_fan!Method
sort_fan!(fan, fan_triangles, tri::Triangulation)

Given a set of triangles in a fan, fan_triangles, associated with some triangulation tri, places all the triangle vertices and sorts them counter-clockwise, placing the results into fan.

source
DelaunayTriangulation.split_marked_vertices!Method
split_marked_vertices!(fan_triangles, tri::Triangulation, marked_vertices)

Given a set of marked_vertices indicating a crossed triangle (like in Figure 9 of this paper), finds all triangles whose three vertices are all in marked_vertices and places them into fan_triangles.

source
DelaunayTriangulation.split_segment!Method
split_segment!(tri::Triangulation, segment, collinear_segments)
split_segment!(segments, segment, collinear_segments)

Splits segment at the segments in collinear_segments, which are assumed to be collinear with segment.

Arguments

  • tri::Triangulation: The Triangulation.
  • segments: The underlying set of segments. This is get_interior_segments(tri) if tri is a Triangulation.
  • segment: The segment to split.
  • collinear_segments: The segments that are collinear with segment.

Outputs

There is no output, as segments is updated in-place.

Example

julia> using DelaunayTriangulation

julia> segments = Set(((2, 3), (3, 5), (10, 12)))
Set{Tuple{Int64, Int64}} with 3 elements:
  (2, 3)
  (3, 5)
  (10, 12)

julia> collinear_segments = [(2, 10), (11, 15), (2, 3)]
3-element Vector{Tuple{Int64, Int64}}:
 (2, 10)
 (11, 15)
 (2, 3)

julia> segment = (3, 5)
(3, 5)

julia> DelaunayTriangulation.split_segment!(segments, segment, collinear_segments)
Set{Tuple{Int64, Int64}} with 4 elements:
  (2, 10)
  (2, 3)
  (11, 15)
  (10, 12)
source
DelaunayTriangulation.triangulate_cavity_cdt!Method
triangulate_cavity_cdt!(tri::Triangulation, V, marked_vertices; rng::Random.AbstractRNG=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel())

Triangulates the cavity V left behind when deleting triangles intersected in a triangulation by an edge, updating tri to do so.

Arguments

  • tri::Triangulation: The Triangulation to update. This should be an empty triangulation.
  • V: The list of polygon vertices, given as a counter-clockwise list of vertices, defining the cavity.
  • tri_fan::Triangulation: The Triangulation to use for the fan of triangles to be re-triangulated. This should be an empty triangulation.
  • marked_vertices: Cache for marking vertices to re-triangulate during the triangulation.
  • fan_triangles: A cache used for sorting and identifying triangles in a fan for retriangulation.

Keyword Arguments

  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use or setup_cavity_cdt.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, but tri is updated in-place.

source
DelaunayTriangulation.get_boundary_edge_mapMethod
get_boundary_edge_map(tri::Triangulation, ij) 
get_boundary_edge_map(tri::Triangulation, i, j)

Returns the value from the key (i, j) in the boundary edge map of tri. The returned value is a Tuple (position, index) so that boundary_nodes = get_boundary_nodes(tri, position) are the boundary nodes associated with the section that (i, j) resides on, and i = get_boundary_nodes(boundary_nodes, index) and j = get_boundary_nodes(boundary_nodes, index + 1).

source
DelaunayTriangulation.contains_boundary_edgeMethod
contains_boundary_edge(tri::Triangulation, ij) -> Bool 
contains_boundary_edge(tri::Triangulation, i, j) -> Bool

Returns true if the boundary edge (i, j) is in tri, and false otherwise. Orientation matters here.

source
DelaunayTriangulation.delete_boundary_node!Method
delete_boundary_node!(tri::Triangulation, pos)

Deletes the boundary node at the specified position pos in tri.

Arguments

  • tri::Triangulation: The Triangulation.
  • pos: The position to delete the node at, given as a Tuple so that delete_boundary_node!(tri, pos) is the same as deleteat!(get_boundary_nodes(tri, pos[1]), pos[2]).

Outputs

There are no outputs, but the boundary nodes of tri are updated in-place.

source
DelaunayTriangulation.get_boundary_nodesMethod
get_boundary_nodes(tri, mnℓ...)

Given a triangulation tri, returns the specified component of the boundary nodes. There are several forms for the methods:

  1. get_boundary_nodes(tri, m): If tri has multiple curves, this returns the mth curve. If tri has multiple sections, this returns the mth section. Otherwise, this returns the mth boundary node.
  2. get_boundary_nodes(tri, m, n): If tri has multiple curves, this returns the nth section of the mth curve. Otherwise, if tri has multiple sections, this returns the nth boundary node of the mth section.
  3. get_boundary_nodes(tri, (m, n)): This is equivalent to get_boundary_nodes(tri, m, n).
  4. get_boundary_nodes(tri::A, ::A): This just returns boundary_nodes.
source
DelaunayTriangulation.get_left_boundary_nodeMethod
get_left_boundary_node(tri::Triangulation, k, ghost_vertex) -> Vertex

Returns the boundary node to the left of the vertex k in tri.

Arguments

  • tri::Triangulation: The Triangulation.
  • k: The boundary vertex.
  • ghost_vertex: The ghost vertex associated with the boundary section that k is on.

Outputs

  • : The vertex left of k on the boundary.
source
DelaunayTriangulation.get_right_boundary_nodeMethod
get_right_boundary_node(tri::Triangulation, k, ghost_vertex) -> Vertex

Returns the boundary node to the right of the vertex k in tri.

Arguments

  • tri::Triangulation: The Triangulation.
  • k: The boundary vertex.
  • ghost_vertex: The ghost vertex associated with the boundary section that k is on.

Outputs

  • r: The vertex right of k on the boundary.
source
DelaunayTriangulation.insert_boundary_node!Method
insert_boundary_node!(tri::Triangulation, pos, node)

Inserts a boundary node into tri at the specified position.

Arguments

  • tri::Triangulation: The Triangulation.
  • pos: The position to insert the node at, given as a Tuple so that insert_boundary_node!(tri, pos, node) is the same as insert!(get_boundary_nodes(tri, pos[1]), pos[2], node).
  • node: The node to insert.

Outputs

There are no outputs, but the boundary nodes of tri are updated in-place.

source
DelaunayTriangulation.split_boundary_edge_at_collinear_segments!Method
split_boundary_edge_at_collinear_segments!(tri::Triangulation, collinear_segments)

Splits a boundary edge into pieces defined by collinear_segments. In particular, if r = collinear_segments and

u = initial(r[1])
v = terminal(r[end]),

then the boundary edge is (u, v) and the edges are split so that all segments in collinear_segments appear instead.

source
DelaunayTriangulation.compute_representative_points!Method
compute_representative_points!(tri::Triangulation; use_convex_hull=!has_boundary_nodes(tri), precision=one(number_type(tri)))

Computes a new set of representative points for tri.

Arguments

  • tri::Triangulation: The Triangulation for which to compute the representative points.

Keyword Arguments

  • use_convex_hull=!has_boundary_nodes(tri): If true, then the representative points are computed using the convex hull of the triangulation. Otherwise, the representative points are computed using the boundary nodes of the triangulation.
  • precision=one(number_type(tri)): The precision to use when computing the representative points via pole_of_inaccessibility.

Output

There are no outputs as tri is updated in-place, but for each curve the representative point is computed using pole_of_inaccessibility.

Exterior curves

While get_exterior_curve_indices(tri) does store the curves corresponding to exterior curves, this function still treats the first curve as the most important exterior curve, computing the representative point so that it is in no holes. In particular, other exterior curves might have representative points that are in a hole of one of their interior holes. This isn't much of a problem, indeed it wouldn't be a significant problem even if we had the representative point in a hole of the first curve, but it is something to be aware of.

source
DelaunayTriangulation.contains_segmentMethod
contains_segment(tri::Triangulation, ij) -> Bool 
contains_segment(tri::Triangulation, i, j) -> Bool

Returns true if (i, j) is a segment in tri, and false otherwise. Both (i, j) and (j, i) are checked.

source

Weighted Triangulations

DelaunayTriangulation.get_distance_to_witness_planeMethod

" getdistancetowitnessplane([kernel::AbstractPredicateKernel = AdaptiveKernel(), ] tri::Triangulation, i, V; cache = nothing) -> Number

Computes the distance between the lifted companion of the vertex i and the witness plane to the triangle V. If V is a ghost triangle and i is not on its solid edge, then the distance is -Inf if it is below the ghost triangle's witness plane and Inf if it is above. If V is a ghost triangle and i is on its solid edge, then the distance returned is the distance associated with the solid triangle adjoining V.

In general, the distance is positive if the lifted vertex is above the witness plane, negative if it is below, and zero if it is on the plane.

The kernel argument determines how this result is computed, and should be one of ExactKernel, FastKernel, and AdaptiveKernel (the default). See the documentation for more information about these choices.

The cache keyword argument is passed to [point_position_relative_to_circumcircle]. Please see the documentation for that function for more information.

See also point_position_relative_to_witness_plane and get_distance_to_plane.

source
DelaunayTriangulation.get_lifted_pointMethod
get_lifted_point(tri::Triangulation, i) -> NTuple{3, Number}

Returns the lifted companion of the ith vertex of tri, in particular (x, y, x^2 + y^2 - w), where w is the ith weight of tri and (x, y) is the ith point of tri.

source
DelaunayTriangulation.get_power_distanceMethod
get_power_distance(tri::Triangulation, i, j) -> Number

Returns the power distance between vertices i and j, defined by ||pᵢ - pⱼ||^2 - wᵢ - wⱼ, where wᵢ and wⱼ are the respective weights.

source
DelaunayTriangulation.get_weightMethod
get_weight(weights, i) -> Number

Gets the ith weight from weights. The default definition for this is weights[i], but this can be extended - e.g., ZeroWeight uses get_weight(weights, i) = 0.0.

If i is not an integer, then the default definition is is_point3(i) ? i[3] : zero(number_type(weights)).

source
DelaunayTriangulation.get_weighted_nearest_neighbourFunction
get_weighted_nearest_neighbour(tri::Triangulation, i, j = rand(each_solid_vertex(tri))) -> Vertex

Using a greedy search, finds the closest vertex in tri to the vertex i (which might not already be in tri), measuring distance in lifted space (i.e., using the power distance - see get_power_distance). The search starts from the vertex j which should be in tri.

source
DelaunayTriangulation.is_submergedFunction
is_submerged([kernel::AbstractPredicateKernel = AdaptiveKernel(), ] tri::Triangulation, i; cache = nothing) -> Bool 
is_submerged([kernel::AbstractPredicateKernel = AdaptiveKernel(), ] tri::Triangulation, i, V; cache = nothing) -> Bool

Returns true if the vertex i is submerged in tri and false otherwise. In the second method, V is a triangle containing tri.

The kernel argument determines how this result is computed, and should be one of ExactKernel, FastKernel, and AdaptiveKernel (the default). See the documentation for more information about these choices.

The cache keyword argument is passed to [point_position_relative_to_circumcircle]. Please see the documentation for that function for more information.

source
DelaunayTriangulation.is_weightedMethod
is_weighted(weights) -> Bool

Returns true if weights represents a set of weights that are not all zero, and false otherwise. Note that even for vectors like zeros(n), this will return true; by default, false is returned only for weights = ZeroWeight().

source

Mesh Refinement

Here are some functions involved with mesh refinement.

DelaunayTriangulation._split_subsegment_curve_bounded!Method
_split_subsegment_curve_bounded!(tri::Triangulation, args::RefinementArguments, e)

Splits a subsegment e of tri at a position determined by split_subcurve! for curve-bounded domains. See split_subsegment!. See also _split_subsegment_curve_bounded_standard! and _split_subsegment_curve_bounded_small_angle!, as well as the original functions _split_subsegment_curve_standard! and _split_subcurve_complex!, respectively, used during boundary enrichment.

source
DelaunayTriangulation.assess_triangle_qualityMethod
assess_triangle_quality(tri::Triangulation, args::RefinementArguments, T) -> Float64, Bool

Assesses the quality of a triangle T of tri according to the RefinementArguments.

Arguments

Output

  • ρ: The radius-edge ratio of the triangle.
  • flag: Whether the triangle is bad quality.

A triangle is bad quality if it does not meet the area constraints, violates the custom constraint, or if it is skinny but neither seditious or nestled.

source
DelaunayTriangulation.check_for_invisible_steiner_pointMethod
check_for_invisible_steiner_point(tri::Triangulation, V, T, flag, c) -> Point, Triangle

Determines if the Steiner point c's insertion will not affect the quality of T, and if so instead changes c to be T's centroid.

Arguments

  • tri::Triangulation: The Triangulation to split a triangle of.
  • V: The triangle that the Steiner point is in.
  • T: The triangle that the Steiner point is from.
  • flag: A Certificate which is Cert.On if the Steiner point is on the boundary of V, Cert.Outside if the Steiner point is outside of V, and Cert.Inside if the Steiner point is inside of V.
  • c: The Steiner point.

Output

  • c′: The Steiner point to use instead of c, which is T's centroid if c is not suitable.
  • V′: The triangle that the Steiner point is in, which is T if c is not suitable.
source
DelaunayTriangulation.check_for_steiner_point_on_segmentMethod
check_for_steiner_point_on_segment(tri::Triangulation, V, V′, new_point, flag, predicates::AbstractPredicateKernel) -> Bool

Checks if the Steiner point with vertex new_point is on a segment. If so, then its vertex is pushed into the offcenter-split list from args, indicating that it should no longer be regarded as a free vertex (see is_free).

Arguments

  • tri::Triangulation: The Triangulation.
  • V: The triangle that the Steiner point was originally in prior to check_for_invisible_steiner_point.
  • V′: The triangle that the Steiner point is in.
  • new_point: The vertex associated with the Steiner point.
  • flag: A Certificate which is Cert.On if the Steiner point is on the boundary of V, Cert.Outside if the Steiner point is outside of V, and Cert.Inside if the Steiner point is inside of V.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

  • onflag: Whether the Steiner point is on a segment or not.
source
DelaunayTriangulation.check_steiner_point_precisionMethod
check_steiner_point_precision(tri::Triangulation, T, c) -> Bool

Checks if the Steiner point c of a triangle T of tri can be computed without precision issues, returning true if there are precision issues and false otherwise.

source
DelaunayTriangulation.compute_split_positionMethod
compute_split_position(tri::Triangulation, args::RefinementArguments, e) -> NTuple{2, Float}

Computes the position to split a segment e of tri at in split_subsegment!.

Arguments

Output

  • mx, my: The position to split the segment at.

This point is computed according to a set of rules:

  1. If e is not a subsegment, meaning it is an input segment, then its midpoint is returned.
  2. If e is a subsegment and the segment adjoins two other distinct segments (one for each vertex) at an acute angle, as determined by segment_vertices_adjoin_other_segments_at_acute_angle, then the point is returned so that e can be split such that one of the new subsegments has a power-of-two length between 1/4 and 1/2 of the length of e, computed using compute_concentric_shell_quarternary_split_position.
  3. If e is a subsegment and the segment adjoins one other segment at an acute angle, as determined by segment_vertices_adjoin_other_segments_at_acute_angle, then the point is returned so that e can be split such that one of the new subsegments has a power-of-two length between 1/3 and 2/3 of the length of e, computed using compute_concentric_shell_ternary_split_position.
  4. Otherwise, the midpoint is returned.
source
DelaunayTriangulation.enqueue_newly_encroached_segments!Method
enqueue_newly_encroached_segments!(args::RefinementArguments, tri::Triangulation) -> Bool

Enqueues all segments that are newly encroached upon after a point insertion into the triangulation into args.queue.

Arguments

Output

  • any_encroached: Whether any segments were newly encroached upon.
source
DelaunayTriangulation.finalise!Method
finalise!(tri::Triangulation, args::RefinementArguments)

Finalises the triangulation after refinement, e.g. by deleting ghost triangles and unlocking the convex hull if needed.

source
DelaunayTriangulation.get_steiner_pointMethod
get_steiner_point(tri::Triangulation, args::RefinementArguments, T) -> Certificate, Point

Computes the Steiner point for a triangle T of tri to improve its quality in split_triangle!.

Arguments

Output

  • precision_flag: A Certificate which is Cert.PrecisionFailure if the Steiner point could not be computed due to precision issues, and Cert.None otherwise.
  • c: The Steiner point. If is_precision_failure(precision_flag), then this is just an arbitrary point of T to ensure type stability.
source
DelaunayTriangulation.is_triangle_nestledMethod
is_triangle_nestled(tri::Triangulation, T, idx) -> Bool

Determines if a triangle T of tri is nestled in the corner of a small input angle.

See also is_triangle_seditious.

Arguments

  • tri::Triangulation: The Triangulation.
  • T: The triangle.
  • idx: The index of the smallest edge of the triangle, so that 1 means uv is the smallest edge, 2 means vw is the smallest edge, and 3 means wu is the smallest edge.

Output

  • flag: Whether the triangle is nestled in the corner of a small input angle.

A triangle is nestled in the corner of a small input angle if it is nestled in the corner of a small input angle and the shortest edge is seditious. The size of the angle is not checked by this function, and is instead determined by assess_triangle_quality.

source
DelaunayTriangulation.is_triangle_seditiousMethod
is_triangle_seditious(tri::Triangulation, args, u, v, w, smallest_idx) -> Bool

Determines if a triangle uvw of tri is seditious according to the RefinementArguments.

See also is_triangle_nestled.

Arguments

  • tri::Triangulation: The Triangulation.
  • args: The RefinementArguments.
  • u, v, w: The vertices of the triangle.
  • smallest_idx: The index of the smallest edge of the triangle, so that 1 means uv is the smallest edge, 2 means vw is the smallest edge, and 3 means wu is the smallest edge.

Output

  • flag: Whether the triangle is seditious.

A triangle is seditious if it is nestled in the corner of a small input angle, or if it is nestled in the corner of a small input angle and the shortest edge is seditious. Here, 'small' is defined by args.constraints.seditious_angle.

source
DelaunayTriangulation.locate_steiner_pointMethod
locate_steiner_point(tri::Triangulation, args::RefinementArguments, T, c) -> Triangle, Cert

Locates the Steiner point c of a triangle T of tri in get_steiner_point. The Steiner point is located by walking from the initial vertex init to c using find_triangle.

Arguments

  • tri::Triangulation: The Triangulation to split a triangle of.
  • args::RefinementArguments: The RefinementArguments for the refinement.
  • T: The triangle that the Steiner point is from.
  • c: The Steiner point.

Output

  • V: The triangle that the Steiner point is in.
  • flag: A Certificate which is Cert.On if the Steiner point is on the boundary of V, Cert.Outside if the Steiner point is outside of V, and Cert.Inside if the Steiner point is inside of V.
source
DelaunayTriangulation.refine!Method
refine!(tri::Triangulation; kwargs...)

Refines the given Triangulation tri to meet the given quality constraints.

See the documentation for more information about mesh refinement, e.g. convergence issues and issues with small input-angles.

Arguments

Keyword Arguments

  • min_angle=30.0: The minimum angle constraint, in degrees.
  • max_angle=180.0: The maximum angle constraint, in degrees.
Maximum angle constraints

Maximum angle constraints are not currently implemented.

  • min_area=get_area(tri) / 1e9: The minimum area constraint.
  • max_area=typemax(number_type(tri)): The maximum area constraint.
  • max_points=max(1_000, num_solid_vertices(tri))^2: The maximum number of vertices allowed in the triangulation. Note that this refers to num_solid_vertices, not the amount returned by num_points.
  • seditious_angle=20.0: The angle at which a triangle is considered seditious, in degrees. See is_triangle_seditious.
  • custom_constraint=(tri, T) -> false: A custom constraint function that takes a Triangulation and a triangle, and returns true if the triangle should be refined and false otherwise.
  • use_circumcenter=true: Whether to insert circumcenters for refining a triangle or generalised Steiner points.
Generalised Steiner points

Generalised Steiner points are not yet implemented. Thus, this argument must be true (and the steiner_scale keyword below is ignored).

  • use_lens=true: Whether to use the diametral lens or the diametral circle for checking encroachment.
  • steiner_scale=0.999: The perturbation factor to use for generalised Steiner points if use_circumcenter=false. (Not currently used - see above.)
  • rng=Random.default_rng(): The random number generator to use in case it is needed during point location.
  • concavity_protection=false: Whether to use concavity protection or not for find_triangle. Most likely not needed, but may help in pathological cases.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

The triangulation is refined in-place.

Duplicate points and unused points

During refinement, points are often deleted, which may often lead to points in get_points(tri) that do not appear anywhere in the triangulation. (This is why we recommend e.g. each_solid_vertex over each_point.) Similarly, since points are deleted, when two triangles have a common circumcenter it might happen (if they are near an input segment) that a point is duplicated inside get_points(tri), in case one circumcenter was deleted previously.

source
DelaunayTriangulation.segment_vertices_adjoin_other_segments_at_acute_angleFunction
segment_vertices_adjoin_other_segments_at_acute_angle(tri::Triangulation, e) -> Int, Vertex

Determines if the vertices of a segment e of tri adjoin other segments at an acute angle.

Arguments

Output

  • num_adjoin: The number of vertices of e that adjoin other segments at an acute angle.
  • adjoin_vert: The vertex of e that adjoins another segment at an acute angle if num_adjoin == 1, and otherwise.
source
DelaunayTriangulation.split_triangle!Method
split_triangle!(tri::Triangulation, args::RefinementArguments, T) -> Certificate

Splits a bad triangle T of tri to improve its quality.

Arguments

Output

  • cert: A Certificate indicating whether the split was successful or not. In particular, returns one of:
    • Cert.SuccessfulInsertion: The triangle was split successfully.
    • Cert.EncroachmentFailure: The triangle was not split successfully as the newly inserted point encroached upon a segment.
    • Cert.PrecisionFailure: The triangle was not split successfully due to precision issues.
source

RTree

The RTrees we work with during boundary enrichment have several functions associated with them.

DelaunayTriangulation.EnlargementValuesType
EnlargementValues

Type for representing the values used in the minimisation of enlargement.

Fields

  • enlargement::Float64: The enlargement of the bounding box of the child being compared with.
  • idx::Int: The index of the child being compared with.
  • area::Float64: The area of the child being compared with.
  • bounding_box::BoundingBox: The bounding box being compared with for enlargement.

Constructor

EnlargementValues(enlargement, idx, area, bounding_box)
EnlargementValues(bounding_box)
source
Base.append!Method
append!(node::AbstractNode, child)

Appends child to node's children. Also updates node's bounding box.

source
Base.delete!Method
delete!(tree::RTree, id_bounding_box::DiametralBoundingBox)

Deletes id_bounding_box from tree.

source
Base.insert!Method
insert!(node::AbstractNode, child, tree::RTree) -> Bool

Inserts child into node in tree. Returns true if the tree's bounding boxes had to be adjusted and false otherwise.

source
Base.insert!Method
insert!(tree::RTree, bounding_box[, level = 1]) -> Bool

Inserts bounding_box into tree. Returns true if the tree's bounding boxes had to be adjusted and false otherwise.

source
Base.iterateMethod
iterate(itr::RTreeIntersectionIterator, state...)

Iterate over the next state of itr to find more intersections with the bounding box in RTreeIntersectionIterator.

source
DelaunayTriangulation.collapse_after_deletion!Method
collapse_after_deletion!(node::AbstractNode, tree::RTree, detached)

Condenses tree after a deletion of one of node's children. The detached argument will contain the nodes that were detached from tree during the condensing process.

source
DelaunayTriangulation.detach!Method
detach!(node::AbstractNode, idx) -> Bool

Detaches the idxth child of node. Returns true if the node's bounding box had to be adjusted and false otherwise.

source
DelaunayTriangulation.find_bounding_boxMethod
find_bounding_box(tree::RTree, id_bounding_box::DiametralBoundingBox) -> Tuple{Leaf{Branch}, Int}

Returns the leaf node and the index in the leaf node's children that id_bounding_box is associated with.

source
DelaunayTriangulation.get_intersectionsMethod
get_intersections(tree::RTree, point::NTuple{2,<:Number}; cache_id=1) -> RTreeIntersectionIterator

Returns an RTreeIntersectionIterator over the elements in tree that intersect with point, representing point as a BoundingBox with zero width and height centered at point. cache_id must be 1 or 2, and determines what cache to use for the intersection query.

source
DelaunayTriangulation.get_next_childMethod
get_next_child(node::AbstractNode, start_idx, need_tests, itr::RTreeIntersectionIterator) -> Int, QueryResult

Returns the index of the next child of node that intersects with the bounding box in itr and the QueryResult of the intersection.

source
DelaunayTriangulation.minimise_enlargementMethod
minimise_enlargement(node, bounding_box) -> EnlargementValues

Returns the EnlargementValues associated with the child of node that minimises the enlargement, where enlargement is defined as the difference between the area of bounding_box and the area of the child's bounding box.

source
DelaunayTriangulation.replace!Method
replace!(node::AbstractNode, left, right, original_bounding_box, tree::RTree)

Replaces the node in tree with left and right. Returns true if the tree's bounding boxes had to be adjusted and false otherwise.

source
DelaunayTriangulation.split!Method
split!(node::AbstractNode, tree::RTree) -> Branch, Branch

Splits node into two other nodes using a linear splitting rule. Returns the two new nodes.

source
DelaunayTriangulation.update_bounding_box!Method
update_bounding_box!(node, idx, original_bounding_box, tree)

Updates the bounding box of node to be the union of its children's bounding boxes. If node has a parent, updates the bounding box of the parent as well if needed.

source

Point Location

Here are some functions related to point location.

DelaunayTriangulation.brute_force_searchMethod
brute_force_search(tri::Triangulation, q; itr = each_triangle(tri), predicates::AbstractPredicateKernel=AdaptiveKernel())

Searches for the triangle containing the point q by brute force. An exception will be raised if no triangle contains the point.

See also find_triangle.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The point to be located.

Keyword Arguments

  • itr = each_triangle(tri): The iterator over the triangles of the triangulation.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

  • V: The triangle containing the point q.
source
DelaunayTriangulation.brute_force_search_enclosing_circumcircleFunction
brute_force_search_enclosing_circumcircle(tri::Triangulation, i, predicates::AbstractPredicateKernel=AdaptiveKernel(); cache = nothing) -> Triangle

Searches for a triangle in tri containing the vertex i in its circumcircle using brute force. If tri is a weighted Delaunay triangulation, the triangle returned instead has the lifted vertex i below its witness plane. If no such triangle exists, (0, 0, 0) is returned. You can control the method used for computing predicates via the predicates argument.

The cache argument is passed to [point_position_relative_to_circumcircle] and should be one of

  • nothing: No cache is used.
  • get_incircle_cache(tri): The cache stored inside tri.
  • AdaptivePredicates.incircleadapt_cache(number_type(tri)): Compute a new cache.

The cache is only needed if an AdaptiveKernel() is used.

source
DelaunayTriangulation.check_for_intersections_with_adjacent_boundary_edgesMethod
check_for_intersections_with_adjacent_boundary_edges(tri::Triangulation, k, q, ghost_vertex=𝒢, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Certificate, Certificate, Vertex, Certificate, Certificate)

Given a boundary vertex k, find a triangle adjacent to k to locate a triangle or edge containing q.

See also search_down_adjacent_boundary_edges, which uses this function to determine an initial direction to search along a straight boundary in case q is collinear with it.

Arguments

  • tri::Triangulation: The Triangulation.
  • k: The boundary vertex to start from.
  • q: The query point.
  • ghost_vertex=𝒢: The ghost vertex corresponding to the boundary that k resides on.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • direction_cert: The direction of q relative to the vertex k along the boundary, given as a Certificate Left, Right, or Outside. If is_outside(direction_cert), then q is not collinear with either of the adjacent boundary edges.
  • q_pos_cert: The position of q relative to the vertex k along the boundary, given as a Certificate Left, Right, On, Outside, or Degenerate. This is similar to direction_cert in that it will be Outside whenever direction_cert is, but this certificate can also be On to indicate that not only is q in the direction given by direction_cert, but it is directly on the edge in that direction. If is_degnerate(q_pos_cert), then q = get_point(tri, next_vertex).
  • next_vertex: The next vertex along the boundary in the direction of q, or k if q is not collinear with either of the adjacent boundary edges.
  • right_cert: The Certificate for the position of q relative to the boundary edge right of k.
  • left_cert: The Certificate for the position of q relative to the boundary edge left of k.
source
DelaunayTriangulation.check_for_intersections_with_interior_edges_adjacent_to_boundary_vertexMethod
check_for_intersections_with_interior_edges_adjacent_to_boundary_vertex(tri::Triangulation, k, q, right_cert, left_cert, predicates::AbstractPredicateKernel=AdaptiveKernel(), store_history=Val(false), history=nothing, ghost_vertex=𝒢) -> (Bool, Vertex, Vertex, Certificate, Certificate)

Checks for intersections between the line pq, where p = get_point(tri, k), and the edges neighbouring p, assuming k is a boundary node. This function should only be used after using check_for_intersections_with_adjacent_boundary_edges.

Arguments

Outputs

The output takes the form (i, j, edge_cert, triangle_cert). Rather than defining each output individually, here are the possible froms of the output:

  • (i, j, Single, Outside): The line pq intersects the edge pᵢpⱼ and (j, i, k) is a positively oriented triangle so that pᵢ is left of pq and pⱼ is right of pq.
  • (i, j, None, Inside): The point q is inside the positively oriented triangle (i, j, k).
  • (0, 0, None, Outside): The point q is outside of the triangulation.
  • (i, j, On, Inside): The point q is on the edge pᵢpⱼ, and thus inside the positively oriented triangle (i, j, k).
  • (i, j, Right, Outside):The pointqis collinear with the edgepᵢpⱼ`, but is off of it and further into the triangulation.
Non-convex geometries

This function assumes that the geometry is convex.

Extended help

This function works in two stages. Firstly, using check_for_intersections_with_single_interior_edge_adjacent_to_boundary_vertex, we check for the intersection of pq with the edges neighbouring the vertex k, rotating counter-clockwise until we find an intersection or reach the other side of the boundary, starting from the first edge counter-clockwise away from the boundary edge right of the vertex k. By keeping track of the positions of pq relative to the current vertex and the previous, we can identify when an intersection is found. If no intersection is found before reaching the boundary edge left of k, then check_for_intersections_with_triangle_left_to_boundary_vertex is used to check the remaining triangle.

source
DelaunayTriangulation.compare_distanceMethod
compare_distance(current_dist, current_idx, pts, i, qx, qy) -> (Number, Vertex)

Computes the minimum of the distance between the ith point of pts and (qx, qy) and current_dist.

Arguments

  • current_dist: The current value for the distance to the point (qx, qy).
  • current_idx: The point of pts corresponding to the distance current_dist.
  • pts: The point set.
  • i: The vertex to compare with current_idx.
  • qx: The x-coordinate of the query point.
  • qy: The y-coordinate of the query point.

Outputs

  • current_dist: The minimum of the distance between the ith point of pts and (qx, qy) and current_dist.
  • current_idx: The point of pts corresponding to the distance current_dist, which will be either i or current_idx.
source
DelaunayTriangulation.concavity_protection_checkFunction
concavity_protection_check(tri::Triangulation, concavity_protection, V, q, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Bool

Check whether the find_triangle algorithm needs to restart. This is only needed if tri is not convex.

Arguments

  • tri::Triangulation: The Triangulation.
  • concavity_protection: Whether this check is needed.
  • V: The triangle that the algorithm has found.
  • q: The query point.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • need_to_restart: Whether the algorithm needs to restart. Will also be false if concavity_protection.

Extended help

This function uses dist to determine whether the query point q is inside or outside of the polygon defined by the triangulation, and also checks the position of q relative to V via point_position_relative_to_triangle. If q is outside of this triangle, then need_to_restart = true. If q is inside this triangle, then issues can still arise due to overlapping ghost triangles from the non-convexity. Thus, depending on the result from dist and whether V is a ghost triangle, need_to_restart will be set accordingly.

source
DelaunayTriangulation.default_num_samplesMethod
default_num_samples(n) -> Integer

Given a number of points n, returns ∛n rounded up to the nearest integer. This is the default number of samples used in the jump-and-march algorithm.

source
DelaunayTriangulation.exterior_find_triangleFunction
exterior_find_triangle(tri::Triangulation, k, q, ghost_vertex=𝒢, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Vertex, Vertex)

Given a query point q outside of the triangulation, find the ghost triangle containing q.

Arguments

  • tri: The Triangulation.
  • k: The exterior ghost vertex to start from.
  • q: The query point.
  • ghost_vertex=𝒢: The ghost vertex corresponding to the boundary that k resides on.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • i: The first vertex of the edge on the boundary adjoining the positively oriented ghost triangle.
  • j: The second vertex of the edge on the boundary adjoining the positively oriented ghost triangle.
Non-convex geometries and internal query points

This function assumes that the geometry is convex. If the geometry is not convex, the returned value may not be correct and should be checked separately. Additionally, if q is actually inside the triangulation, then the result is meaningless.

Extended help

This function works by first finding the position of q relative to pₘp, where pₘ is the representative point for the ghost_vertex and p = get_point(tri, k). Depending on this position, we rotate counter-clockwise if q is left of the line using exterior_find_triangle_rotate_left and clockwise otherwise using exterior_find_triangle_rotate_right. By keeping track of the current position of q and its position relative to the next ghost edge, we can identify when q resides inside a ghost triangle.

source
DelaunayTriangulation.find_triangleMethod
find_triangle(tri, q; kwargs...) -> Triangle[, Bool]

Find the triangle in the triangulation tri containing the query point q using the jump-and-march algorithm.

Ghost triangles

For this function to work best, the triangulation should have ghost triangles, which you can add using add_ghost_triangles! in case tri does not already have them. Without ghost triangles, the function may not be able to find the correct triangle containing q.

For the variables defined below, you may want to refer to the extended help which also gives some warnings and notes.

Arguments

Keyword Arguments

  • predicates::AbstractPredicateKernel=ExactKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • point_indices=each_solid_vertex(tri): The indices of the vertices to consider as possible starting points for the algorithm.
  • m=default_num_samples(num_vertices(point_indices)): The number of samples to use when selecting the initial point.
  • try_points=(): A list of points to try as the initial point in addition to the m sampled.
  • rng=Random.default_rng(): The random number generator to use.
  • k=select_initial_point(tri, q; point_indices, m, try_points, rng): The initial point to start the algorithm from. See select_initial_point.
  • store_history=Val(false): Whether to store the history of the algorithm.
  • history=nothing: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • maxiters=2 + num_exterior_curves(tri) - num_solid_vertices(tri) + num_solid_edges(tri): The maximum number of iterations to perform before restarting the algorithm with restart_find_triangle.
  • concavity_protection=false: Whether to use concavity protection. See concavity_protection_check. This is only needed if your triangulation is not convex.
  • use_barriers::Val{U}=Val(false): Whether to stop searching beyond any segments in the triangulation.

Output

  • V: The triangle containing q, with type given by triangle_type(tri).

If you have use_barriers == Val(true), then we also return

  • invisible_flag: false if the triangle was found without hitting a barrier, and true otherwise.

Extended help

The algorithm underlying this function is complicated and broken into many parts. Here, we describe a brief overview of the algorithm, but note that the documentation contains a much more detailed description.

  1. Firstly, the algorithm is initialised depending on whether k is a boundary or an interior vertex, using initialise_find_triangle_boundary_vertex or initialise_find_triangle_interior_vertex respectively.
  2. From the initial triangle (i, j, k) chosen, we then check if q is one of pᵢ, pⱼ, and p = pₖ and then return according to find_triangle_return_on_vertex if needed.
  3. If we do not return above, we need to step from the initial triangle towards q. Since we put pᵢ and pⱼ to the left and right of the line pq, respectively, this means that we step until the triangle pᵢpⱼq is no longer positively oriented. So, while the triangle is positively oriented, we step according to find_triangle_across_triangle.
  4. If we have not yet returned and the triangle is no longer positively oriented, we check if the triangle is degenerate using find_triangle_degenerate_arrangement and reinitialise the algorithm if needed. Otherwise, we have found the triangle containing q and return the triangle.

Here are some additional warnings and notes for the variables defined in this function.

Restarting the algorithm

If the algorithm restarts, then the initial point k is selected again using select_initial_point, and the algorithm is restarted from there. This is done if the algorithm gets stuck in a loop, or if the algorithm is not able to find the correct triangle containing q after maxiters iterations. For a convex geometry, maxiters can be safely ignored, as the sequence of triangles visited is acyclic [see H. Edelsbrunner, An acyclicity theorem for cell complexes in d dimensions, Combinatorica 10 (1990) 251-260.)].

Found triangles with barriers

If you are using barriers, it will be your responsibility to verify that any found triangle from this function actually contains the triangle. This can be verified using the returned flag (see below), although the point might still be on the triangle's boundary.

Walking past vertices of barriers

If you are using barriers, it is possible that the algorithm can walk past vertices of barriers. One way this can happen is if the initial search line intersects a vertex, in which case the segment might not be considered. Another way this can happen is if you start the algorithm directly on a segment vertex, in which case the algorithm can go past it (e.g. this means that it is possible that a ghost triangle might still be returned if you start the algorithm on a boundary node).

Some notes about the output:

Hitting barriers

If a barrier is hit before any initial triangle is properly identified, the returned triangle is (0, 0, 0); this is only possible if use_barriers == Val(true). Moreover, if use_barriers == Val(true), the final triangle may not even be valid if invisible_flag == true (defined below).

Non-convex geometries

While this function does still work for non-convex geometries, it may be significantly slower than for convex geometries, as most of the details of the algorithm assume that the geometry is convex, and so the algorithm may have to restart many times at new initial vertices k.

source
DelaunayTriangulation.find_triangle_across_triangleMethod
find_triangle_across_triangle(tri::Triangulation, q, k, predicates, store_history, history, maxiters, cur_iter, concavity_protection, arrangement, original_k, last_changed, p, i, j, pᵢ, pⱼ) -> (Bool, Bool, Bool, Triangle, Integer, Certificate, Vertex, Vertex, Vertex, Point, Point, Integer, Integer)

Walks across the current triangle past the edge (i, j). progressing the find_triangle algorithm.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point.
  • k: The vertex that the algorithm started from.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • store_history: Whether to store the history of the algorithm.
  • history: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • maxiters: The maximum number of iterations to perform before restarting the algorithm with restart_find_triangle.
  • cur_iter: The current iteration of the algorithm.
  • concavity_protection: Whether to use concavity protection. See concavity_protection_check. This is only needed if your triangulation is not convex.
  • arrangement: A Certificate defining the orientation of the triangle pᵢpⱼq.
  • original_k: The original vertex that the algorithm started from.
  • last_changed: The last vertex that was changed in the algorithm.
  • p: The point corresponding to the vertex original_k.
  • i: The first vertex of the triangle adjoining k to step from.
  • j: The second vertex of the triangle adjoining k to step from.
  • pᵢ: The point corresponding to the vertex i.
  • pⱼ: The point corresponding to the vertex j.

Outputs

  • restart_flag: Whether the algorithm needs to be restarted.
  • return_flag: Whether the algorithm can return immediately, returning V.
  • reinitialise_flag: Whether the algorithm needs to be reinitialised at a new vertex k. (This would only be needed if !has_ghost_triangles(tri).)
  • V: The triangle stepped into.
  • cur_iter: The new number of iterations of the algorithm.
  • arrangement: A Certificate defining the orientation of the triangle pᵢpⱼq with the updated values of i and j.
  • k: The new vertex that the algorithm is at.
  • last_changed: The last vertex that was changed in the algorithm.
  • original_k: The original vertex that the algorithm started from.
  • pᵢ: The point corresponding to the vertex i.
  • pⱼ: The point corresponding to the vertex j.
  • i: The first vertex of the triangle adjoining k to step from.
  • j: The second vertex of the triangle adjoining k to step from.

Extended help

This part of the algorithm is relatively complicated because there are many cases that need to be accounted for. Here we give a brief description of how this step works, and note that the documentation contains a much more detailed description.

  1. Firstly, we need to check whether k is an exterior ghost vertex or not. If k is an exterior ghost vertex, then this means that we are stepping outside of the triangulation. Thus, we use exterior_find_triangle to find where q is, starting from the last_changed vertex. If concavity_protection = true, then concavity_protection_check is used to determine if a restart is needed, or if we can return safely. If we reach this step but !has_ghost_triangles(tri), then the algorithm should need to be reinitialised since q should not be outside of the triangulation, and so we return with reinitialise_flag = true.
  2. Now we consider the case where k is not an exterior ghost vertex. We move forward by updating the value of k so that k = get_adjacent(tri, i, j), and then consider where pₖ is relative to the line pq.

2a. If pₖ is to the right of pq, then we should update j by j = k, ensuring that j is always to the right of pq.

2b. If pₖ is to the left of pq, then we should update i by i = k, ensuring that i is always to the left of pq.

2c. The alternative to 2a and 2b is that pₖ is collinear with the edge of pq, which could mean that q is in the current triangle or it is in a triangle further away. We compute a Certificate that determines where q is relative to the triangle pᵢpⱼpₖ. If q is inside or on this triangle, then we return, restarting if necessary according to concavity_protection and concavity_protection_check. If we do not yet need to return, then we need to make a decision as to which of i and j to update, noting that we want i to be left of pq and j to be right of pq, but this is no longer unambiguous since pₖ is collinear with pq. We make this decision according to last_changed: If last_changed = i, then moving left is what caused us to find this collinear edge, and so we send k left by letting i = k. Otherwise, we send k right by letting j = k.

  1. Now having stepped forward, we recompute the Certificate for arrangement and return, setting restart_flag = true if cur_iters ≥ maxiters.
source
DelaunayTriangulation.find_triangle_degenerate_arrangementMethod
find_triangle_degenerate_arrangement(tri::Triangulation, q, k, predicates::AbstractPredicateKernel, store_history::F, history, pᵢ, pⱼ, i, j) -> Bool

Given a degenerate arrangement pᵢpⱼq, reinitialise the jump and march algorithm.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point.
  • k: The vertex that the algorithm started from.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • store_history: Whether to store the history of the algorithm.
  • history: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • pᵢ: The point corresponding to the vertex i.
  • pⱼ: The point corresponding to the vertex j.
  • i: The first vertex of the triangle adjoining k to step from.
  • j: The second vertex of the triangle adjoining k to step from.

Outputs

  • Bool: Whether the algorithm needs to be restarted.
source
DelaunayTriangulation.find_triangle_return_on_vertexMethod
find_triangle_return_on_vertex(tri::Triangulation, q, k, p, pᵢ, pⱼ, i, j, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Bool, Bool, Triangle)

Check if q is one of the vertices of the triangle (i, j, k) and return if needed.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point.
  • k: The vertex k that the algorithm started from.
  • p: The point corresponding to the vertex k.
  • pᵢ: The point corresponding to the vertex i.
  • pⱼ: The point corresponding to the vertex j.
  • i: The first vertex of the triangle adjoining k to start from.
  • j: The second vertex of the triangle adjoining k to start from.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • restart_flag: Whether the algorithm needs to be restarted.
  • return_flag: Whether the algorithm can return immediately, returning V.
  • V: The triangle (i, j, k).

Extended help

An extra check is made in this algorithm for the case that the point that q is equal to is one of the points corresponding to a ghost vertex, so it may be for example that q == pᵢ but is_ghost_vertex(i), in which case the algorithm would need to restart.

source
DelaunayTriangulation.initialise_find_triangle_boundary_vertexMethod
initialise_find_triangle_boundary_vertex(tri::Triangulation, q, k, predicates::AbstractPredicateKernel, store_history, history, ghost_vertex, concavity_protection) -> (Bool, Bool, Triangle, Point, Vertex, Vertex, Point, Point)

Initialise the jump-and-march algorithm for a boundary vertex k.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point.
  • k: The boundary vertex to start from.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • store_history: Whether to store the history of the algorithm.
  • history: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • ghost_vertex: The ghost vertex corresponding to the boundary that k resides on.
  • concavity_protection: Whether to use concavity protection. See concavity_protection_check. This is only needed if your triangulation is not convex.

Outputs

  • restart_flag: Whether the algorithm needs to be restarted.
  • return_flag: Whether the algorithm can return immediately, returning V.
  • V: Either a triangle that can returned if return_flag = true, or some triangle that is used for type stability for this return value.
  • p: The point corresponding to the vertex k, or it may be q if the algorithm is going to be restarted or return_flag = true.
  • i: The first vertex of the triangle adjoining k to start from, or k if the algorithm is going to be restarted or return_flag = true.
  • j: The second vertex of the triangle adjoining k to start from, or k if the algorithm is going to be restarted or return_flag = true.
  • pᵢ: The point corresponding to the vertex i, or it may be q if the algorithm is going to be restarted or return_flag = true.
  • pⱼ: The point corresponding to the vertex j, or it may be q if the algorithm is going to be restarted or return_flag = true.

Extended help

There are multiple stages to this initialisation, starting from check_for_intersections_with_adjacent_boundary_edges.

  • If it is found that q is not outside of the triangulation, so that q is collinear with one of the boundary edges, then we use search_down_adjacent_boundary_edges to find where to start, noting that we can return immediately if q is found to be on an adjacent boundary edge. Otherwise, exterior_find_triangle can then be used to find the ghost triangle containing q; if concavity_protection = true, then concavity_protection_check is used to determine if a restart is needed.

  • If is is not found that q is outside of the triangulation yet based on information from the adjacent boundary edges, then we need to check the neighbouring interior edges using check_for_intersections_with_interior_edges_adjacent_to_boundary_vertex, returning early if q is found to be inside one of the neighbouring triangles. If the line pq, where p = get_point(tri, k), does not intersect any of the neighbouring edges and it is not inside any of the neighbouring triangles, then it must be outside of the triangulation and so we use exterior_find_triangle to find the triangle; as before, concavity_protection_check is used on the found ghost triangle if needed. If there is an intersection, then we return the triangle containing the intersection point that we can start the algorithm from, and its associated vertices and points.

source
DelaunayTriangulation.initialise_find_triangle_interior_vertexMethod
initialise_find_triangle_interior_vertex(tri::Triangulation, q, k, predicates::AbstractPredicateKernel, store_history::F, history, rng) -> (Bool, Point, Vertex, Vertex, Point, Point)

Initialise the jump-and-march algorithm for an interior vertex k.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point.
  • k: The interior vertex to start from.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • store_history: Whether to store the history of the algorithm.
  • history: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • rng: The random number generator to use.

Outputs

  • restart_flag: Whether the algorithm needs to be restarted.
  • p: The point corresponding to the vertex k.
  • i: The first vertex of the triangle adjoining k to start from.
  • j: The second vertex of the triangle adjoining k to start from.
  • pᵢ: The point corresponding to the vertex i.
  • pⱼ: The point corresponding to the vertex j.

Extended help

This function works by simply using select_initial_triangle_interior_vertex to find the initial triangle to start from. A check is made to see if the edge (i, j) refers to a non-existent edge (0, 0), in which case the algorithm needs to be restarted.

source
DelaunayTriangulation.prepare_initial_edgeFunction
prepare_initial_edge(tri::Triangulation, edges, p, q, rng::Random.AbstractRNG=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Vertex, Vertex, Point, Point, Certificate, Certificate)

Selects a random edge from the set of edges edges and computes the certificates for the points corresponding to the initial and terminal vertices of the edge.

Arguments

  • tri::Triangulation: The Triangulation.
  • edges: The set of edges to sample from.
  • p: The initial point.
  • q: The query point.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • i: The initial vertex of the edge.
  • j: The terminal vertex of the edge.
  • pᵢ: The point corresponding to i.
  • pⱼ: The point corresponding to j.
  • line_cert_i: The Certificate for pᵢ's position relative to the oriented line pq.
  • line_cert_j: The Certificate for pⱼ's position relative to the oriented line pq.
source
DelaunayTriangulation.restart_find_triangleMethod
restart_find_triangle(tri::Triangulation, q, store_history, history, rng, maxiters, cur_iter, concavity_protection, num_restarts, use_barriers, predicates::AbstractPredicateKernel) -> Triangle[, Bool]

Restart the find_triangle algorithm, or use brute_force_search to find q if num_restarts ≥ 25.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point.
  • store_history: Whether to store the history of the algorithm.
  • history: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • rng: The random number generator to use.
  • maxiters: The maximum number of iterations to perform before restarting the algorithm with restart_find_triangle.
  • cur_iter: The current iteration of the algorithm.
  • concavity_protection: Whether to use concavity protection. See concavity_protection_check. This is only needed if your triangulation is not convex.
  • num_restarts: The number of times the algorithm has been restarted.
  • use_barriers: Whether to use barriers, stopping the algorithm at any segment.
  • predicates::AbstractPredicateKernel: Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • V: The triangle containing q.

In addition, if use_barriers = Val(true), then a second output is returned, which is a boolean indicating whether the algorithm reached a barrier (true) or not (false).

source
DelaunayTriangulation.search_down_adjacent_boundary_edgesMethod
search_down_adjacent_boundary_edges(tri::Triangulation, k, q, direction_cert, q_pos_cert, next_vertex, predicates::AbstractPredicateKernel=AdaptiveKernel(), store_history=Val(false), history=nothing, ghost_vertex=𝒢) -> (Bool, Certificate, Vertex, Vertex, Vertex)

Starting at the boundary vertex k, walks down the boundary in the direction of q until finding q or finding that it is outside of the triangulation.

See also check_for_intersections_with_adjacent_boundary_edges, which uses this function to determine an initial direction to search along.

Arguments

Outputs

  • return_flag: Whether to return, or throw an exception.
  • q_pos_cert: A Certificate that is On if q is on the edge (u, v), and Outside if q is outside of the triangulation.
  • u: If is_on(q_pos_cert), this is the first vertex of a positively oriented triangle that q is on, so that q is on the edge (u, v). Otherwise, (u, v, w) is a ghost triangle close to q.
  • v: If is_on(q_pos_cert), this is the second vertex of a positively oriented triangle that q is on, so that q is on the edge (u, v). Otherwise, (u, v, w) is a ghost triangle close to q.
  • w: If is_on(q_pos_cert), this is the third vertex of a positively oriented triangle that q is on, so that q is on the edge (u, v) and w = get_adjacent(tri, u, v). Otherwise, (u, v, w) is a ghost triangle close to q.
Non-convex geometries

This function assumes that the geometry is convex. The function will still be able to return, but is_outside(q_pos_cert) may not necessarily mean q is outside of the triangulation. The main function find_triangle will have to restart the algorithm if it is found that is_outside(q_pos_cert) was incorrect.

Extended help

This function works by stepping along vertices on the boundaries in the direction specified by direction_cert, using search_right_down_adjacent_boundary_edges if is_right(direction_cert) and search_left_down_adjacent_boundary_edges otherwise. In these functions, a while loop is used to keep stepping until q_pos_cert, which is updated at each iteration, changes value.

source
DelaunayTriangulation.select_initial_pointFunction
select_initial_point(tri::Triangulation, q; kwargs...) -> Vertex

Selects the initial point for find_triangle to start from.

Arguments

  • tri::Triangulation: The Triangulation.
  • q: The query point. Can be either a point or a vertex - if it is a vertex, the corresponding point get_point(tri, q) will be used.

Keyword Arguments

  • point_indices=each_solid_vertex(tri): The indices to sample from.
  • m=default_num_samples(num_vertices(point_indices)): The number of samples to take. Replacement is not used, so there may be duplicates.
  • try_points=(): A list of points to try in addition to those randomly sampled.
  • allow_boundary_points=!is_disjoint(tri): Whether to allow boundary points to be selected.
  • rng=Random.default_rng(): The random number generator to use.

Outputs

  • i: The index of the point closest to q out of those queried.
source
DelaunayTriangulation.select_initial_triangle_interior_vertexMethod
select_initial_triangle_interior_vertex(tri::Triangulation, 
    k, 
    q, 
    predicates::AbstractPredicateKernel=AdaptiveKernel(),
    store_history=Val(false), 
    history=nothing, 
    rng::Random.AbstractRNG=Random.default_rng()) -> (Point, Vertex, Vertex, Point, Point)

Selects the initial triangle for find_triangle to start from, for the case where k is an interior vertex.

Arguments

  • tri::Triangulation: The Triangulation.
  • k: The vertex to start from.
  • q: The query point.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • store_history=Val(false): Whether to store the history of the algorithm.
  • history=nothing: The history of the algorithm. If store_history, then this should be a PointLocationHistory object.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • p: The point corresponding to k.
  • i: The initial vertex of the triangle.
  • j: The terminal vertex of the triangle.
  • pᵢ: The point corresponding to i.
  • pⱼ: The point corresponding to j.
Non-convex geometries

This function assumes that the geometry is convex. To deal with this, when an infinite loop is detected we return for both i and j, and then let find_triangle handle how to correct the algorithm from there.

Extended help

This part of the algorithm works by rotating around the vertex k, looking for a triangle whose edges adjoining k are to the left and to the right of k. By choosing the initial edge at random via prepare_initial_edge, and computing the position of q relative to this initial edge, the rotation will be either clockwise or counter-clockwise, and the triangle is then found using either select_initial_triangle_clockwise or select_initial_triangle_counterclockwise, respectively.

In case the initial edge is collinear with the line pq, where p = get_point(tri, q), then fix_initial_collinear_edge_for_interior_vertex to find a non-collinear edge resample more edges from prepare_initial_edge if necessary.

source
DelaunayTriangulation.select_random_edgeFunction
select_random_edge(tri::Triangulation, edges, rng::Random.AbstractRNG=Random.default_rng()) -> (Vertex, Vertex, Point, Point)

Selects a random edge from the set of edges edges.

Arguments

  • tri::Triangulation: The Triangulation.
  • edges: The set of edges to sample from.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use.

Outputs

  • i: The initial vertex of the edge.
  • j: The terminal vertex of the edge.
  • pᵢ: The point corresponding to i.
  • pⱼ: The point corresponding to j.
source

Triangulation Operations

Some of the triangulation operations have internal functions associated with them.

DelaunayTriangulation.delete_holes!Method
delete_holes!(tri::Triangulation)

Deletes all the exterior faces to the boundary nodes specified in the triangulation tri.

Extended help

This function works in several stages:

  1. First, find_all_points_to_delete is used to identify all points in the exterior faces.
  2. Once all the points to delete have been found, all the associated triangles are found using find_all_triangles_to_delete, taking care for any incorrectly identified triangles and points.
  3. Once the correct set of triangles to delete has been found, they are deleted using delete_all_exterior_triangles!.
source
DelaunayTriangulation.find_all_points_to_delete!Method
find_all_points_to_delete!(points_to_process, tri::Triangulation, seed, all_bn)

Starting at seed, finds more points to spread to and mark for deletion.

Arguments

  • points_to_process: The current list of points marked for deletion.
  • tri::Triangulation: The Triangulation.
  • seed: The seed vertex to start spreading from.
  • all_bn: All the boundary nodes in the triangulation, obtained from get_all_boundary_nodes.

Outputs

There are no outputs, as points_to_process is updated in-place.

Extended help

This function works by considering the neighbours around the vertex seed. For each neighbouring vertex, we designate that as a new seed, and consider if it needs to be added into points_to_process according to its distance from the triangulation computed from distance_to_polygon. We then call find_all_points_to_delete! recursively again on the new seed.

source
DelaunayTriangulation.find_all_points_to_deleteMethod
find_all_points_to_delete(tri::Triangulation) -> Set{Int}

Returns a set of all the points that are in the exterior faces of the triangulation tri.

Extended help

This function works by 'spreading' from some initial vertex. In particular, starting at each boundary node, we spread outwards towards adjacent vertices, recursively spreading so that all exterior points are identified with the help of find_all_points_to_delete!.

source
DelaunayTriangulation.find_all_triangles_to_deleteMethod
find_all_triangles_to_delete(tri::Triangulation, points_to_process) -> Set{Triangle}

Returns a set of all the triangles that are in the exterior faces of the triangulation tri.

Arguments

Outputs

  • triangles_to_delete: The set of triangles that are in the exterior faces of the triangulation tri.

Extended help

This function works in two stages.

  1. Firstly, all the non-boundary vertices, i.e. those from points_to_process, are processed. For each vertex v, the triangles adjoining it, given by get_adjacent2vertex(tri, v), aremarked for deletion.
  2. Next, all the boundary vertices need to be processed and carefully analysed to determine if any other triangles need to be deleted since, for example, a triangle may be adjoining three vertices that are all boundary vertices, and it might not be obvious if it is inside or outside of the triangulation. By applying dist to compute the distance between the triangle's centroid and the triangulation, the triangle can be accurately marked for deletion if it is outside of the triangulation.
source
DelaunayTriangulation.check_delete_point_argsMethod
check_delete_point_args(tri::Triangulation, vertex, S) -> Bool

Checks that the vertex vertex can be deleted from the triangulation tri. Returns true if so, and throws an InvalidVertexDeletionError otherwise. This will occur if:

  1. vertex is a boundary node of tri.
  2. vertex is a ghost vertex of tri.
  3. vertex adjoins a segment of tri.
source
DelaunayTriangulation.delete_point!Method
delete_point!(tri::Triangulation, vertex; kwargs...)

Deletes the vertex of tri, retriangulating the cavity formed by the surrounding polygon of vertex using triangulate_convex.

It is not possible to delete vertices that are on the boundary, are ghost vertices, or adjoin a segment of tri. See also check_delete_point_args.

Point deletion

This function will not actually delete the corresponding coordinates from get_points(tri), nor will it remove the associated weight from get_weights(tri).

Arguments

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • store_event_history=Val(false): Whether to store the event history of the triangulation from deleting the point.
  • event_history=nothing: The event history of the triangulation from deleting the point. Only updated if store_event_history is true, in which case it needs to be an InsertionEventHistory object.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator to use for the triangulation.
source
DelaunayTriangulation.get_surrounding_polygonMethod
get_surrounding_polygon(tri::Triangulation, u; skip_ghost_vertices=false) -> Vector

Returns the counter-clockwise sequence of neighbours of u in tri.

Arguments

Keyword Arguments

  • skip_ghost_vertices=false: Whether to skip ghost vertices in the returned polygon.

Outputs

  • S: The surrounding polygon. This will not be circular, meaning S[begin] ≠ S[end]. In case u is an exterior ghost vertex, the returned polygon is a clockwise list of vertices for the associated boundary curve. If you do not have ghost triangles and you try to get the surrounding polygon of a ghost vertex, then this function may return an invalid polygon.
source
DelaunayTriangulation.complete_split_edge_and_legalise!Function
complete_split_edge_and_legalise!(tri::Triangulation, i, j, r, store_event_history=Val(false), event_history=nothing; predicates::AbstractPredicateKernel=AdaptiveKernel())

Given a triangulation tri, an edge (i, j), and a point r, splits both (i, j) and (j, i) at r using split_edge! and then subsequently legalises the new edges with legalise_split_edge!.

Arguments

  • tri::Triangulation: The Triangulation.
  • i: The first vertex of the edge to split.
  • j: The second vertex of the edge to split.
  • r: The vertex to split the edge at.
  • store_event_history=Val(false): Whether to store the event history of the flip.
  • event_history=nothing: The event history. Only updated if store_event_history is true, in which case it needs to be an InsertionEventHistory object.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, as tri is updated in-place.

source
DelaunayTriangulation.legalise_split_edge!Function
legalise_split_edge!(tri::Triangulation, i, j, k, r, store_event_history=Val(false), event_history=nothing; predicates::AbstractPredicateKernel=AdaptiveKernel())

Legalises the newly added edges in tri after the edge (i, j) was split using split_edge!.

See also complete_split_edge_and_legalise!.

Arguments

  • tri::Triangulation: The Triangulation.
  • i: The first vertex of the edge that was split.
  • j: The second vertex of the edge that was split.
  • k: The vertex that was originally adjacent to (i, j).
  • r: The vertex that (i, j) was split at.
  • store_event_history=Val(false): Whether to store the event history of the flip.
  • event_history=nothing: The event history. Only updated if store_event_history is true, in which case it needs to be an InsertionEventHistory object.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, as tri is updated in-place.

source
DelaunayTriangulation.split_edge!Function
split_edge!(tri::Triangulation, i, j, r, store_event_history=Val(false), event_history=nothing)

Splits the edge (i, j) in tri at the vertex r. For the triangulation to be valid after this splitting, it is assumed that r is collinear with, or at least very close to collinear with, the edge (i, j).

See also legalise_split_edge! and complete_split_edge_and_legalise!.

Arguments

  • tri::Triangulation: The Triangulation.
  • i: The first vertex of the edge to split.
  • j: The second vertex of the edge to split.
  • r: The vertex to split the edge at.
  • store_event_history=Val(false): Whether to store the event history of the flip.
  • event_history=nothing: The event history. Only updated if store_event_history is true, in which case it needs to be an InsertionEventHistory object.

Outputs

There is no output, as tri is updated in-place.

Handling unoriented edges

The triangulation will only be updated as if (i, j) has been split rather than also (j, i). You will need to call split_edge! again with (j, i) if you want to split that edge as well.

source
DelaunayTriangulation.complete_split_triangle_and_legalise!Method
complete_split_triangle_and_legalise!(tri::Triangulation, i, j, k, r)

Splits the triangle (i, j, k) at the vertex r, assumed to be inside the triangle, and legalises the newly added edges in tri.

Arguments

  • tri::Triangulation: The Triangulation.
  • i: The first vertex of the triangle.
  • j: The second vertex of the triangle.
  • k: The third vertex of the triangle.
  • r: The vertex to split the triangle at.

Outputs

There is no output, as tri is updated in-place.

source
DelaunayTriangulation.legalise_split_triangle!Method
legalise_split_triangle!(tri::Triangulation, i, j, k, r; predicates::AbstractPredicateKernel=AdaptiveKernel())

Legalises the newly added edges in tri after the triangle (i, j, k) was split using split_triangle!.

See also complete_split_triangle_and_legalise!.

Arguments

  • tri::Triangulation: The Triangulation.
  • i: The first vertex of the triangle.
  • j: The second vertex of the triangle.
  • k: The third vertex of the triangle.
  • r: The vertex to split the triangle at.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There is no output, as tri is updated in-place.

source
DelaunayTriangulation.split_triangle!Method
split_triangle!(tri::Triangulation, i, j, k, r)

Splits the triangle (i, j, k) at the vertex r, assumed to be inside the triangle.

See also legalise_split_triangle! and complete_split_triangle_and_legalise!.

Arguments

  • tri::Triangulation: The Triangulation.
  • i: The first vertex of the triangle.
  • j: The second vertex of the triangle.
  • k: The third vertex of the triangle.
  • r: The vertex to split the triangle at.

Outputs

There is no output, but tri will be updated so that it now contains the triangles (i, j, r), (j, k, r), and (k, i, r).

source
DelaunayTriangulation._safe_get_adjacentMethod
_safe_get_adjacent(tri::Triangulation, uv) -> Vertex

This is the safe version of get_adjacent, which is used when the triangulation has multiple sections, ensuring that the correct ghost vertex is returned in case uv is a ghost edge.

source
DelaunayTriangulation.get_adjacentMethod
get_adjacent(tri::Triangulation, uv) -> Vertex
get_adjacent(tri::Triangulation, u, v) -> Vertex

Returns the vertex w such that (u, v, w) is a positively oriented triangle in the underlying triangulation, or if no such triangle exists.

source
DelaunayTriangulation.get_neighboursMethod
get_neighbours(tri::Triangulation, u) -> Set{Vertex}

Returns the set of neighbours of u in tri. Note that, if has_ghost_triangles(tri), then some of the neighbours and vertices will be ghost vertices.

source
DelaunayTriangulation.get_neighboursMethod
get_neighbours(tri::Triangulation) -> Dict{Vertex, Set{Vertex}}

Returns the neighbours map of tri. Note that, if has_ghost_triangles(tri), then some of the neighbours and vertices will be ghost vertices.

source
DelaunayTriangulation.num_neighboursMethod
num_neighbours(tri::Triangulation, u) -> Integer

Returns the number of neighbours of u in tri. Note that, if has_ghost_triangles(tri), then some of the neighbours counted might be ghost vertices if u is a boundary vertex.

source
DelaunayTriangulation.get_pointMethod
get_point(tri::Triangulation, i) -> NTuple{2, Number}
get_point(tri::Triangulation, i...) -> NTuple{length(i), NTuple{2, Number}}

Returns the coordinates corresponding to the vertices i... of tri, given as a Tuple of the form (x, y) for each point. If i is a ghost vertex, then the coordinates of the representative point of the curve associated with i are returned instead.

source
DelaunayTriangulation.num_pointsMethod
num_points(tri::Triangulation) -> Integer

Returns the number of points in tri.

Danger

If tri has vertices that are not yet present in the triangulation, e.g. if you have deleted vertices or have some submerged vertices in a weighted triangulation, then the corresponding points will still be counted in this function. It is recommended that you instead consider num_vertices, num_solid_vertices, or num_ghost_vertices.

source
DelaunayTriangulation.construct_positively_oriented_triangleFunction
construct_triangle(tri::Triangulation, i, j, k, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Triangle

Returns a triangle in tri from the vertices i, j, and k such that the triangle is positively oriented.

You can use the predicates argument to determine how the orientation predicate is computed. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for more discussion on these choices.

source
DelaunayTriangulation.contains_triangleMethod
contains_triangle(tri::Triangulation, T) -> (Triangle, Bool)
contains_triangle(tri::Triangulation, i, j, k) -> (Triangle, Bool)

Tests whether tri contains T = (i, j, k) up to rotation, returning

  • V: The rotated form of T that is in tri, or simply T if T is not in tri.
  • flag: true if T is in tri, and false otherwise.
source
DelaunayTriangulation.num_trianglesMethod
num_triangles(tri::Triangulation) -> Integer

Returns the number of triangles in tri. Note that, if has_ghost_triangles(tri), then some of these triangles will be ghost triangles.

source

Voronoi Tessellations

Here are some functions related to the computation of unbounded Voronoi tessellations.

DelaunayTriangulation.add_edge_to_voronoi_polygon!Method
add_edge_to_voronoi_polygon!(B, vorn::VoronoiTessellation, i, k, S, m, encountered_duplicate_circumcenter) -> (Vertex, Bool, Vertex)

Add the next edge to the Voronoi polygon for the point i in the VoronoiTessellation vorn.

Arguments

  • B: The vector of circumcenters defining the polygon.
  • vorn: The VoronoiTessellation.
  • i: The polygon index.
  • k: The vertex to add.
  • S: The surrounding polygon of i. See get_surrounding_polygon.
  • m: The index of the next vertex in S.
  • encountered_duplicate_circumcenter: Whether or not a duplicate circumcenter has been encountered.

Outputs

  • ci: The index for the circumcenter of the triangle considered.
  • encountered_duplicate_circumcenter: Whether or not a duplicate circumcenter has been encountered.
  • k: The next vertex in S after the input k.
source
DelaunayTriangulation.add_voronoi_polygon!Method
add_voronoi_polygon!(vorn::VoronoiTessellation, i) -> Vector

Add the Voronoi polygon for the point i to the VoronoiTessellation vorn.

Arguments

Outputs

  • B: The vector of circumcenters defining the polygon. This is a circular vector, i.e. B[begin] == B[end].
source
DelaunayTriangulation.close_voronoi_polygon!Method
close_voronoi_polygon!(vorn::VoronoiTessellation, B, i, encountered_duplicate_circumcenter, prev_ci)

Close the Voronoi polygon for the point i in the VoronoiTessellation vorn.

Arguments

  • vorn: The VoronoiTessellation.
  • B: The vector of circumcenters defining the polygon.
  • i: The polygon index.
  • encountered_duplicate_circumcenter: Whether or not a duplicate circumcenter has been encountered.
  • prev_ci: The previous circumcenter index.

Outputs

There are no outputs, as vorn and B are modified in-place.

source
DelaunayTriangulation.get_next_triangle_for_voronoi_polygonMethod
get_next_triangle_for_voronoi_polygon(vorn::VoronoiTessellation, i, k, S, m) -> (Vertex, Vertex)

Get the next triangle for the Voronoi polygon for the point i in the VoronoiTessellation.

Arguments

Outputs

  • ci: The index for the circumcenter of the next triangle.
  • k: The next vertex in S after the input k.
source

Clipped Voronoi Tessellations

Here are some functions related to the computation of clipped Voronoi tessellations.

DelaunayTriangulation.add_segment_intersection!Method
add_segment_intersection!(segment_intersections, boundary_sites, intersection_point, incident_polygon::I) where {I} -> Integer

Adds the intersection_point into the list of segment_intersections.

Arguments

  • segment_intersections: The list of segment intersections.
  • boundary_sites: A mapping from boundary sites to the indices of the segment intersections that are incident to the boundary site.
  • intersection_point: The intersection point to add.
  • incident_polygon: The index of the polygon that is incident to the intersection point.

Outputs

  • idx: The index of the intersection point in the list of segment intersections. If the intersection point already exists in the list, then the index of the existing point is returned and used instead.
source
DelaunayTriangulation.add_to_intersected_edge_cache!Method
add_to_intersected_edge_cache!(intersected_edge_cache, u, v, a, b)

Add the edge uv to the list of intersected edges.

Arguments

  • intersected_edge_cache: The list of intersected edges.
  • u: The first vertex of the edge of the Voronoi polygon intersecting the edge ab of the boundary.
  • v: The second vertex of the edge of the Voronoi polygon intersecting the edge ab of the boundary.
  • a: The first vertex of the edge of the boundary.
  • b: The second vertex of the edge of the boundary.

Outputs

There are no outputs, as intersected_edge_cache is modified in-place.

source
DelaunayTriangulation.classify_intersections!Method
classify_intersections!(intersected_edge_cache, left_edge_intersectors, right_edge_intersectors, current_edge_intersectors, left_edge, right_edge, current_edge)

Classify the intersections in intersected_edge_cache into left_edge_intersectors, right_edge_intersectors, and current_edge_intersectors based on whether they intersect left_edge, right_edge, or current_edge, respectively.

Arguments

  • intersected_edge_cache: The list of intersected edges currently being considered.
  • left_edge_intersectors: The set of sites that intersect the edge to the left of an edge currently being considered.
  • right_edge_intersectors: The set of sites that intersect the edge to the right of an edge currently being considered.
  • current_edge_intersectors: The set of sites that intersect the current edge being considered.
  • left_edge: The edge to the left of e on the boundary.
  • right_edge: The edge to the right of e on the boundary.
  • current_edge: The edge on the boundary being considered.

Outputs

There are no outputs, but left_edge_intersectors, right_edge_intersectors, or current_edge_intersectors are updated all in-place depending on the type of intersection for each edge in intersected_edge_cache.

source
DelaunayTriangulation.clip_all_polygons!Method
clip_all_polygons!(vorn::VoronoiTessellation, n, boundary_sites, exterior_circumcenters, equal_circumcenter_mapping)

Clip all of the polygons in the Voronoi tessellation.

Arguments

  • vorn: The VoronoiTessellation.
  • n: The number of vertices in the tessellation before clipping.
  • boundary_sites: A dictionary of boundary sites.
  • exterior_circumcenters: Any exterior circumcenters to be filtered out.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.

Outputs

There are no outputs, but the polygons are clipped in-place.

source
DelaunayTriangulation.clip_polygon!Method
clip_polygon!(vorn::VoronoiTessellation, n, points, polygon, new_verts, exterior_circumcenters, equal_circumcenter_mapping)

Clip the polygon polygon by removing the vertices that are outside of the domain and adding the new vertices new_verts to the polygon.

Arguments

  • vorn: The VoronoiTessellation.
  • n: The number of vertices in the tessellation before clipping.
  • points: The polygon points of the tessellation.
  • polygon: The index of the polygon to be clipped.
  • new_verts: The indices of the new vertices that are added to the polygon.
  • exterior_circumcenters: Any exterior circumcenters to be filtered out.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.

Outputs

There are no outputs, but the polygon is clipped in-place.

source
DelaunayTriangulation.clip_voronoi_tessellation!Method
clip_voronoi_tessellation!(vorn::VoronoiTessellation; rng=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel())

Clip the Voronoi tessellation vorn to the convex hull of the generators in vorn.

Arguments

Keyword Arguments

  • clip_polygon=nothing: If clip=true, then this is the polygon to clip the Voronoi tessellation to. If nothing, the convex hull of the triangulation is used. The polygon should be defined as a Tuple of the form (points, boundary_nodes) where the boundary_nodes are vertices mapping to coordinates in points, adhering to the usual conventions for defining boundaries. Must be a convex polygon.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There are no outputs, but the Voronoi tessellation is clipped in-place.

source
DelaunayTriangulation.dequeue_and_process!Function
dequeue_and_process!(vorn, polygon_edge_queue, edges_to_process,
    intersected_edge_cache, left_edge_intersectors, right_edge_intersectors, current_edge_intersectors,
    processed_pairs, boundary_sites, segment_intersections, exterior_circumcenters, equal_circumcenter_mapping,
    rng::Random.AbstractRNG=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel())

Dequeue an edge from polygon_edge_queue and process it. If polygon_edge_queue is empty, then we process the first edge in edges_to_process.

Arguments

  • vorn: The VoronoiTessellation.
  • polygon_edge_queue: The queue of edges that need to be processed.
  • edges_to_process: The edges that need to be processed.
  • intersected_edge_cache: A cache of intersected edges.
  • left_edge_intersectors: The intersection points of left_edge with other edges.
  • right_edge_intersectors: The intersection points of right_edge with other edges.
  • current_edge_intersectors: The intersection points of current_edge with other edges.
  • processed_pairs: A set of pairs of edges and polygons that have already been processed.
  • boundary_sites: A dictionary of boundary sites.
  • segment_intersections: A dictionary of segment intersections.
  • exterior_circumcenters: A dictionary of exterior circumcenters.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.
  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There are no outputs. Instead, the caches and queues are updated in-place.

Extended help

This function works as follows:

  1. Firstly, if there are no edges queued in polygon_edge_queue, then we enqueue the first in edge in edges_to_process using enqueue_new_edge!.
  2. We then dequeue the next edge to be processed. If the edge has already been processed, then we return early.
  3. If we're still here, then we process the (edge, polygon) pair enqueued from polygon_edge_queue using process_polygon!. This function checks for intersections of the edge with the polygon.
  4. Once the polygon has been processed, we then needed to classify all of the intersections using classify_intersections!, which determines, for each intersection, if the intersection is with edge, or with the edge left of edge, or to the edge right of edge.
  5. Then, process_intersection_points! is used to process the intersection points, enqueueing new edges when needed.
  6. We then delete the edge from edges_to_process if it is in there and return.
source
DelaunayTriangulation.enqueue_new_edge!Function
enqueue_new_edge!(polygon_edge_queue, vorn::VoronoiTessellation, e)

Enqueue the edge e of the boundary to be processed.

Arguments

Outputs

There are no outputs, as polygon_edge_queue is modified in-place.

source
DelaunayTriangulation.find_all_exterior_circumcentersMethod
find_all_exterior_circumcenters(vorn::VoronoiTessellation, clip_points, clip_vertices) -> Set{I}

Finds all the polygon vertices in vorn that are outside of the polygon defined by (clip_points, clip_vertices). The return is the set of all exterior vertices.

source
DelaunayTriangulation.find_all_intersecting_polygonsMethod
find_all_intersecting_polygons(vorn::VoronoiTessellation, clip_points, clip_vertices, exterior_circumcenters) -> (Set{I}, Set{I})

Classifies all polygons in vorn into three sets depending on how they intersect the polygon defined by (clip_points, clip_vertices), returning (interior, intersecting, exterior), where

  • interior are the polygons that are entirely inside the polygon defined by (clip_points, clip_vertices),
  • possibly_intersecting are the polygons that might intersect the polygon defined by (clip_points, clip_vertices).

The exterior_circumcenters are the circumcenters that are outside of the domain, see find_all_exterior_circumcenters.

source
DelaunayTriangulation.find_all_intersectionsMethod
find_all_intersections(vorn::VoronoiTessellation; rng=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Dict, Vector, Set, Dict)

Find all intersections between the edges of the Voronoi tessellation and the boundary of the polygon.

Arguments

Keyword Arguments

  • rng::Random.AbstractRNG=Random.default_rng(): The random number generator.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • boundary_sites: A dictionary of boundary sites.
  • segment_intersections: The intersection points.
  • exterior_circumcenters: The circumcenters that are outside of the domain.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.

Extended help

This algorithm works as follows:

  1. First, using initialise_clipping_arrays, we initialise the arrays that we will use to store the intersections, and queue up all boundary edges for processing.
  2. Then, starting with the first edge in edges_to_process, we dequeue an edge from polygon_edge_queue and process it via dequeue_and_process!.
  3. We repeat step 2 until polygon_edge_queue and edges_to_process are both empty.
  4. We then return.
source
DelaunayTriangulation.fix_no_intersections!Method
fix_no_intersections!(vorn::VoronoiTessellation; predicates::AbstractPredicateKernel=AdaptiveKernel())

In the case that the Voronoi tiles have no intersections at all with the convex hull, this function adds in the missing intersections. clip_polygon is used for this.

source
DelaunayTriangulation.get_neighbouring_boundary_edgesMethod
get_neighbouring_boundary_edges(tri::Triangulation, e) -> (Edge, Edge)

Returns the two boundary edges adjacent to the boundary edge e in the triangulation tri.

Arguments

  • tri::Triangulation: a triangulation
  • e: The boundary edge.

Outputs

  • left_e: The left edge.
  • right_e: The right edge.
source
DelaunayTriangulation.get_shared_vertexMethod
get_shared_vertex(e, f) -> Vertex

Returns the vertex shared by the edges e and f, or if they do not share a vertex.

Arguments

  • e: The first edge.
  • f: The second edge.

Outputs

  • u: The shared vertex.

Example

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.get_shared_vertex((1, 3), (5, 7))
0

julia> DelaunayTriangulation.get_shared_vertex((1, 3), (3, 7))
3

julia> DelaunayTriangulation.get_shared_vertex((10, 3), (10, 5))
10

julia> DelaunayTriangulation.get_shared_vertex((9, 4), (9, 5))
9
source
DelaunayTriangulation.initialise_clipping_arraysMethod
initialise_clipping_arrays(vorn::VoronoiTessellation) -> (Set{E}, Queue{Tuple{E,I}}, Dict{I,Set{I}}, NTuple{2,F}[], Set{Tuple{E,I}}, Pair{E,E}[], Set{I}, Set{E}, Set{E}, Set{E}, Dict{I,I})

Initialise the arrays used in the clipping algorithm for the VoronoiTessellation vorn.

Arguments

Outputs

  • edges_to_process: The set of edges that are to be processed.
  • polygon_edge_queue: The queue of edges that are to be processed.
  • boundary_sites: A mapping from boundary sites to the indices of the segment intersections that are incident to the boundary site.
  • segment_intersections: The list of segment intersections.
  • processed_pairs: The set of pairs of edges and polygons that have been processed.
  • intersected_edge_cache: The list of intersected edges currently being considered.
  • exterior_circumcenters: The list of circumcenters of sites that are outside the boundary.
  • left_edge_intersectors: The set of sites that intersect the edge to the left of an edge currently being considered.
  • right_edge_intersectors: The set of sites that intersect the edge to the right of an edge currently being considered.
  • current_edge_intersectors: The set of sites that intersect the current edge being considered.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.
source
DelaunayTriangulation.process_intersection_points!Method
process_intersection_points!(polygon_edge_queue, vorn, current_incident_polygon,
    left_edge_intersectors, right_edge_intersectors, current_edge_intersectors,
    left_edge, right_edge, current_edge, processed_pairs, segment_intersections, boundary_sites)

Process the intersection points in left_edge_intersectors, right_edge_intersectors, and current_edge_intersectors and add the new edges to polygon_edge_queue if necessary. Special care is taken to not miss any corner points.

Arguments

  • polygon_edge_queue: The queue of edges that need to be processed.
  • vorn: The VoronoiTessellation.
  • current_incident_polygon: The index of the current polygon being processed.
  • left_edge_intersectors: The intersection points of left_edge with other edges.
  • right_edge_intersectors: The intersection points of right_edge with other edges.
  • current_edge_intersectors: The intersection points of current_edge with other edges.
  • left_edge: The left edge of the current polygon.
  • right_edge: The right edge of the current polygon.
  • current_edge: The current edge of the current polygon.
  • processed_pairs: A set of pairs of edges and polygons that have already been processed.
  • segment_intersections: A dictionary of segment intersections.
  • boundary_sites: A dictionary of boundary sites.

Outputs

There are no outputs, but the caches and queues are updated in-place.

Extended help

The rules are based on the paper "Efficient Computation of Clipped Voronoi Diagram for Mesh Generation" by Yan, Wang, Levy, and Liu. Namely, an edge that intersects a boundary edge and one adjacent to it has its shared vertex added to the queue together with the current polygon (current_incident_polygon) being considered, and any intersections have the adjacent polygon added to the queue together with the intersecting edge. (These are not strictly the rules in the paper.)

This function works as follows:

  1. First, assuming that there is more than one triangle in the underlying triangulation of vorn, we need to consider left_edge and right_edge individually.
  2. The procedure for each edge is the same, so here we just describe the left_edge. If there are any intersectors with the left_edge, and neither of (left_edge, current_incident_polygon) or (reverse_edge(left_edge), current_incident_polygon) have already been processed (i.e., in processed_pairs), then we enqueue (left_edge, i) and (left_edge, j) into polygon_edge_queue, where i and j are the vertices of left_edge which correspond to polygons. This will ensure that we can find intersections next to this polygon.
  3. After enqueueing these pairs, we also need to protect against corner points, which we check for by considering current_incident_polygon ∈ all_indices, where all_indices are the vertices of left_edge, right_edge, and current_edge. If this is true, and if the shared vertex of current_edge and left_edge is equal to current_incident_polygon, then we need to add the point generator of current_incident_polygon as an intersection. This need comes from having to worry about corners, i.e. points where the two unbounded polygons meet and go directly left and right of a vertex so that that vertex is not considered an intersection; this point needs to be included.
  4. Once the left_edge and right_edge have been processed as above, we need to then consider all of left_edge, right_edge, and current_edge, and each of the intersections through the respective edge. This step is done regardless of whether there is a single triangle in the underlying triangulation. The procedure for each edge is the same, so let us just describe the current_edge. For each edge uv in the current_edge_intersectors, we need to get the polygon adjacent to that edge. Then, if (current_edge, adjacent_incident_polygon) or (reverse_edge(current_edge), adjacent_incident_polygon) have not been processed, we enqueue (current_edge, adjacent_incident_polygon).
  5. Once the edges have all been processed as above, we return.
source
DelaunayTriangulation.process_polygon!Function
process_polygon!(vorn::VoronoiTessellation, e, incident_polygon, boundary_sites, segment_intersections, intersected_edge_cache, exterior_circumcenters, equal_circumcenter_mapping, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Edge, Edge, Edge)

Processes the polygon incident_polygon for all of its intersections based on the boundary edge e.

Arguments

  • vorn::VoronoiTessellation: The VoronoiTessellation.
  • e: The edge on the boundary being considered.
  • incident_polygon: The index of the polygon being considered.
  • boundary_sites: The mapping from the indices of the sites on the boundary to the indices of the edges on the boundary that they intersect.
  • segment_intersections: The list of segment intersections.
  • intersected_edge_cache: A cache of the edges that have been intersected by the ray from u to v.
  • exterior_circumcenters: A list of the circumcenters of the sites that are outside the convex hull of the sites on the boundary.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • left_edge: The edge to the left of e on the boundary.
  • right_edge: The edge to the right of e on the boundary.
  • e: The edge on the boundary being considered.

In addition to these outputs, the caches are also updated in-place.

Extended help

This function works as follows:

  1. First, for the current edge e, we get the edges left_edge and right_edge that neighbour it via get_neighbouring_boundary_edges.
  2. For each edge of the incident_polygon, we need to process it depending on whether the edge (u, v) is finite, between two ghosts, going out to infinity, or coming in from infinity. If the edge is between two ghosts, we skip the edge. For rays that go out or in to infinity, we use process_ray_intersection! and process_ray_intersection_with_other_edges! to process the intersection of the ray with the boundary edges. The function process_ray_intersection_with_other_edges! is needed since rays going out to infinity may have to go through other boundary edges in order to do so, e.g. at a corner it may be that it crosses two boundary edges. For finite segments, process_segment_intersection! is used to process the intersection. We apply this function with each of e, left_edge, and right_edge to check for all intersections.
  3. The function is done once each of the polygon edges has been considered.
source
DelaunayTriangulation.process_ray_intersection!Function
process_ray_intersection!(
    vorn::VoronoiTessellation,
    u,
    v,
    incident_polygon,
    intersected_edge_cache,
    segment_intersections,
    boundary_sites,
    exterior_circumcenters,
    equal_circumcenter_mapping,
    predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Point

Process the intersection of the Voronoi polygon of the site u with the ray emanating from the circumcenter of the site v.

Arguments

  • vorn: The VoronoiTessellation.
  • u: The index of the site u, given as a ghost vertex for the associated ghost triangle.
  • v: The index of the site v.
  • incident_polygon: The index of the Voronoi polygon of the site u that is incident to the ray emanating from the circumcenter of the site v.
  • intersected_edge_cache: The list of intersected edges currently being considered.
  • segment_intersections: The list of segment intersections.
  • boundary_sites: A mapping from boundary sites to the indices of the segment intersections that are incident to the boundary site.
  • exterior_circumcenters: The list of circumcenters of sites that are outside the boundary.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • p: The coordinates of the intersection.

In addition to the point p, add_segment_intersection! is also updated to incorporate the new intersection point, as is add_to_intersected_edge_cache!.

source
DelaunayTriangulation.process_ray_intersection_with_other_edges!Function
process_ray_intersection_with_other_edges!(vorn::VoronoiTessellation,
    u,
    v,
    e,
    left_edge,
    right_edge,
    r,
    segment_intersections,
    boundary_sites,
    incident_polygon,
    equal_circumcenter_mapping,
    intersected_edge_cache,
    predicates::AbstractPredicateKernel=AdaptiveKernel())

Process the intersection of the ray from the ghost site u to the site v with the edges e, left_edge and right_edge.

Arguments

  • vorn::VoronoiTessellation: The VoronoiTessellation.
  • u: The index of the ghost site.
  • v: The index of the site u is going to.
  • e: The edge on the boundary being considered.
  • left_edge: The edge to the left of e on the boundary.
  • right_edge: The edge to the right of e on the boundary.
  • r: The coordinates of the intersection of the ray from u to v with some edge. If any(isnan, r), then the ray does not intersect any edge and we skip.
  • segment_intersections: The list of segment intersections.
  • boundary_sites: The mapping from the indices of the sites on the boundary to the indices of the edges on the boundary that they intersect.
  • incident_polygon: The index of the polygon that contains the intersection of the ray from u to v with the boundary.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.
  • intersected_edge_cache: A cache of the edges that have been intersected by the ray from u to v.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

There are no outputs, but add_segment_intersection! and add_to_intersected_edge_cache! are used to update the intersection objects.

source
DelaunayTriangulation.process_segment_intersection!Function
process_segment_intersection!(
    vorn::VoronoiTessellation,
    u,
    v,
    e,
    incident_polygon,
    intersected_edge_cache,
    segment_intersections,
    boundary_sites,
    exterior_circumcenters,
    equal_circumcenter_mapping,
    predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Point

Process the intersection of the Voronoi polygon's edge (u, v) with the edge e of the boundary, returning the coordinates of the intersection and updating via add_segment_intersection!.

Arguments

  • vorn: The VoronoiTessellation.
  • u: The index of the site u.
  • v: The index of the site v.
  • e: The edge e of the boundary.
  • incident_polygon: The index of the Voronoi polygon currently being considered.
  • intersected_edge_cache: The list of intersected edges currently being considered.
  • segment_intersections: The list of segment intersections.
  • boundary_sites: A mapping from boundary sites to the indices of the segment intersections that are incident to the boundary site.
  • exterior_circumcenters: The list of circumcenters of sites that are outside the boundary.
  • equal_circumcenter_mapping: A mapping from the indices of the segment intersections that are equal to the circumcenter of a site to the index of the site.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • p: The coordinates of the intersection. If there is no intersection, this is (NaN, NaN).

In addition to the point p, add_segment_intersection! is also updated to incorporate the new intersection point, as is add_to_intersected_edge_cache!.

source
DelaunayTriangulation._get_rayFunction
_get_ray(vorn, i, ghost_vertex, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Point, Point)

Extracts the ray from the ith polygon of vorn corresponding to the ghost_vertex, where ghost_vertex here means that get_polygon(vorn, i)[ghost_vertex] is a ghost vertex.

Arguments

  • vorn: The VoronoiTessellation.
  • i: The index of the polygon.
  • ghost_vertex: The index of the ghost vertex in the polygon.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • p: The first point of the ray.
  • q: A second point of the ray, so that pq gives the direction of the ray (which extends to infinity).
source
DelaunayTriangulation.clip_bounded_polygon_to_bounding_boxMethod
clip_bounded_polygon_to_bounding_box(vorn::VoronoiTessellation, i, bounding_box; predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Vector{NTuple{2,Number}}

Clips the ith polygon of vorn to bounding_box.

See also clip_polygon.

Arguments

  • vorn: The VoronoiTessellation.
  • i: The index of the polygon.
  • bounding_box: The bounding box to clip the polygon to.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • coords: The coordinates of the clipped polygon. This is a circular vector.
source
DelaunayTriangulation.clip_unbounded_polygon_to_bounding_boxMethod
clip_unbounded_polygon_to_bounding_box(vorn::VoronoiTessellation, i, bounding_box; predicates::AbsractPredicateType=AdaptiveKernel()) -> Vector{NTuple{2,Number}}

Clips the ith polygon of vorn to bounding_box. The polygon is assumed to be unbounded. See also clip_polygon.

Use the keyword arguments predicates to determine how predicates are computed. Should be one of ExactKernel, AdaptiveKernel, and FastKernel. See the documentation for more information about these choices.

source
DelaunayTriangulation.get_bounded_polygon_coordinatesMethod
get_bounded_polygon_coordinates(vorn::VoronoiTessellation, i, bounding_box; predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Vector{NTuple{2,Number}}

Returns the coordinates of the ith polygon of vorn, clipped to bounding_box.

Use the keyword arguments predicates to determine how predicates are computed. Should be one of ExactKernel, AdaptiveKernel, and FastKernel. See the documentation for more information about these choices.

source
DelaunayTriangulation.get_clipping_poly_structsMethod
get_clipping_poly_structs(vorn::VoronoiTessellation, i, bounding_box) -> (Polygon, Polygon)

Returns the polygons used for clipping the ith polygon of vorn to bounding_box.

See also clip_polygon.

Arguments

  • vorn: The VoronoiTessellation.
  • i: The index of the polygon.
  • bounding_box: The bounding box to clip the polygon to.

Outputs

  • poly: The polygon to clip.
  • clip_poly: The polygon to clip to.
source
DelaunayTriangulation.get_new_polygon_indicesMethod
get_new_polygon_indices(vorn, vertices) -> (Vector{Int}, Vector{NTuple{2,Float64}}, Tuple{Int, Int})

Returns the new vertices and points of the polygon, as well as the indices of the ghost vertices in the polygon.

Arguments

Outputs

  • new_vertices: The new vertices of the polygon. This is not a circular vector. The vertices corresponding to a ghost vertex will be given by the ghost vertex itself.
  • new_points: The new points of the polygon. This is not a circular vector. The points corresponding to a ghost vertex will be given by by (NaN, NaN).
  • ghost_vertices: The indices of the ghost vertices in new_vertices.
source
DelaunayTriangulation.get_polygon_coordinatesFunction
get_polygon_coordinates(vorn::VoronoiTessellation, i, bounding_box=nothing; predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Vector{NTuple{2,Number}}

Returns the coordinates of the polygon with index i in vorn. If bounding_box is provided, then the polygon is clipped to the bounding box. If the polygon is unbounded, then bounding_box must be provided.

See also get_unbounded_polygon_coordinates and get_bounded_polygon_coordinates.

Arguments

  • vorn: The VoronoiTessellation.
  • i: The index of the polygon.
  • bounding_box=nothing: The bounding box to clip the polygon to. If nothing, then the polygon is not clipped. If the polygon is unbounded, then bounding_box must be provided.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • coords: The coordinates of the polygon. This is a circular vector.
source
DelaunayTriangulation.get_unbounded_polygon_coordinatesMethod
get_unbounded_polygon_coordinates(vorn::VoronoiTessellation, i, bounding_box; predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Vector{NTuple{2,Number}}

Returns the coordinates of the ith polygon of vorn, clipped to bounding_box. The polygon is assumed to be unbounded.

Use the keyword arguments predicates to determine how predicates are computed. Should be one of ExactKernel, AdaptiveKernel, and FastKernel. See the documentation for more information about these choices.

source
DelaunayTriangulation.grow_polygon_outside_of_boxFunction
grow_polygon_outside_of_box(vorn::VoronoiTessellation, i, bounding_box, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Vector{Int}, Vector{NTuple{2,Number}})

Truncates the unbounded edges of the ith polygon of vorn so that the line connecting the truncated unbounded edges is entirely outside of bounding_box.

Arguments

  • vorn: The VoronoiTessellation.
  • i: The index of the polygon. The polygon must be unbounded.
  • bounding_box: The bounding box to clip the polygon to. See also polygon_bounds.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • new_vertices: The new vertices of the polygon. This is not a circular vector.
  • new_points: The new points of the polygon. This is not a circular vector.
source
DelaunayTriangulation.is_first_ghost_vertexMethod
is_first_ghost_vertex(cell, i) -> Bool

Assuming that the circular vector cell is such that ghost vertices only appear next to eachother in cell and there are only two, tests if i is the first ghost vertex in cell.

See also is_last_ghost_vertex.

Arguments

  • cell: The circular vector.
  • i: The index of the vertex in cell.

Outputs

  • flag: true if i is the first ghost vertex in cell, and false otherwise.
source
DelaunayTriangulation.is_last_ghost_vertexMethod
is_last_ghost_vertex(cell, i) -> Bool

Assuming that the circular vector cell is such that ghost vertices only appear next to eachother in cell, tests if i is the last ghost vertex in cell.

See also is_first_ghost_vertex.

Arguments

  • cell: The circular vector.
  • i: The index of the vertex in cell.

Outputs

  • flag: true if i is the last ghost vertex in cell, and false otherwise.
source
DelaunayTriangulation.maximum_distance_to_boxMethod
maximum_distance_to_box(a, b, c, d, p) -> Number

Computes the maximum squared distance from the point p to the box with corners (a, c), (b, c), (b, d), (a, d).

Arguments

  • a: The minimum x-coordinate of the box.
  • b: The maximum x-coordinate of the box.
  • c: The minimum y-coordinate of the box.
  • d: The maximum y-coordinate of the box.
  • p: The point.

Outputs

  • dist: The maximum squared distance from p to the box.
source
DelaunayTriangulation.liang_barskyMethod
liang_barsky(a, b, c, d, p, q) -> (Point, Point)

Applies the Liang-Barsky algorithm to find the intersection of the line segment pq with the rectangle [a, b] × [c, d].

Arguments

  • p: The first point of the line segment.
  • q: The second point of the line segment.
  • a: The minimum x-coordinate of the rectangle.
  • b: The maximum x-coordinate of the rectangle.
  • c: The minimum y-coordinate of the rectangle.
  • d: The maximum y-coordinate of the rectangle.

Output

  • u: The first coordinate of the intersection, or (NaN, NaN) if there is no intersection.
  • v: The second coordinate of the intersection, or (NaN, NaN) if there is no intersection.
source
DelaunayTriangulation.clip_polygonMethod
clip_polygon(vertices, points, clip_vertices, clip_points; predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Vector

Clip a polygon defined by (vertices, points) to a convex clip polygon defined by (clip_vertices, clip_points) with the Sutherland-Hodgman algorithm. The polygons should be defined in counter-clockwise order.

Arguments

  • vertices: The vertices of the polygon to be clipped.
  • points: The underlying point set that the vertices are defined over.
  • clip_vertices: The vertices of the clipping polygon.
  • clip_points: The underlying point set that the clipping vertices are defined over.

Keyword Arguments

  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Output

  • clipped_polygon: The coordinates of the clipped polygon, given in counter-clockwise order and clipped_polygon[begin] == clipped_polygon[end].
source

Centroidal Voronoi Tessellations

Here are some functions related to the computation of centroidal Voronoi tessellations.

DelaunayTriangulation._centroidal_smooth_itrFunction
_centroidal_smooth_itr(vorn::VoronoiTessellation, points, rng, predicates::AbstractPredicateKernel=AdaptiveKernel(); kwargs...) -> (VoronoiTessellation, Number)

Performs a single iteration of the centroidal smoothing algorithm.

Arguments

  • vorn: The VoronoiTessellation.
  • points: The underlying point set. This is a deepcopy of the points of the underlying triangulation.
  • rng: The random number generator.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Keyword Arguments

  • clip_polygon=nothing: If clip=true, then this is the polygon to clip the Voronoi tessellation to. If nothing, the convex hull of the triangulation is used. The polygon should be defined as a Tuple of the form (points, boundary_nodes) where the boundary_nodes are vertices mapping to coordinates in points, adhering to the usual conventions for defining boundaries.
  • kwargs...: Extra keyword arguments passed to retriangulate.

Outputs

source
DelaunayTriangulation.centroidal_smoothMethod
centroidal_smooth(vorn::VoronoiTessellation; maxiters=1000, tol=default_displacement_tolerance(vorn), rng=Random.default_rng(), predicates::AbstractPredicateKernel=AdaptiveKernel(), kwargs...) -> VoronoiTessellation

Smooths vorn into a centroidal tessellation so that the new tessellation is of a set of generators whose associated Voronoi polygon is that polygon's centroid.

Arguments

Keyword Arguments

  • maxiters=1000: The maximum number of iterations.
  • clip_polygon=nothing: If clip=true, then this is the polygon to clip the Voronoi tessellation to. If nothing, the convex hull of the triangulation is used. The polygon should be defined as a Tuple of the form (points, boundary_nodes) where the boundary_nodes are vertices mapping to coordinates in points, adhering to the usual conventions for defining boundaries. Must be a convex polygon.
  • tol=default_displacement_tolerance(vorn): The displacement tolerance. See default_displacement_tolerance for the default.
  • rng=Random.default_rng(): The random number generator.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • kwargs...: Extra keyword arguments passed to retriangulate.

Outputs

Extended help

The algorithm is simple. We iteratively smooth the generators, moving them to the centroid of their associated Voronoi polygon for the current tessellation, continuing until the maximum distance moved of any generator is less than tol. Boundary generators are not moved.

source
DelaunayTriangulation.default_displacement_toleranceMethod
default_displacement_tolerance(vorn::VoronoiTessellation) -> Number

Returns the default displacement tolerance for the centroidal smoothing algorithm. The default is given by max_extent / 1e4, where max_extent = max(width, height), where width and height are the width and height of the bounding box of the underlying point set.

Arguments

Outputs

  • tol: The default displacement tolerance.
source
DelaunayTriangulation.move_generator_to_centroid!Method
move_generator_to_centroid!(points, vorn::VoronoiTessellation, generator) -> Number

Moves the generator generator to the centroid of its Voronoi cell. Returns the distance moved.

Arguments

  • points: The underlying point set. This is a deepcopy of the points of the underlying triangulation.
  • vorn: The VoronoiTessellation.
  • generator: The generator to move.

Outputs

  • dist: The distance moved.
source

Triangulating Curve-Bounded Domains

We have several functions related to the triangulation of curve-bounded domains.

DelaunayTriangulation.coarse_discretisation!Method
coarse_discretisation!(points, boundary_nodes, boundary_curve; n=0)

Constructs an initial coarse discretisation of a curve-bounded domain with bonudary defines by (points, boundary_nodes, boundary_curves), where boundary_nodes and boundary_curves should come from convert_boundary_curves!. The argument n is the amount of times to split an edge. If non-zero, this should be a power of two (otherwise it will be rounded up to the next power of two). If it is zero, then the splitting will continue until the maximum total variation over any subcurve is less than π/2.

source
DelaunayTriangulation.compute_split_positionFunction
compute_split_position(enricher::BoundaryEnricher, i, j, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> (Float64, Float64, NTuple{2,Float64})

Gets the point to split the edge (i, j) at.

Arguments

  • enricher::BoundaryEnricher: The enricher.
  • i: The first point of the edge.
  • j: The second point of the edge.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • t: The parameter value of the split point.
  • Δθ: The total variation of the subcurve (i, t). If a split was created due to a small angle, this will be set to zero.
  • ct: The point to split the edge at.
source
DelaunayTriangulation.enrich_boundary!Method
enrich_boundary!(enricher::BoundaryEnricher; predicates::AbstractPredicateKernel=AdaptiveKernel())

Enriches the initial boundary defined inside enricher, implementing the algorithm of Gosselin and Ollivier-Gooch (2007). At the termination of the algorithm, all edges will contain no other points inside their diametral circles.

The predicates argument determines how predicates are computed, and should be one of ExactKernel, FastKernel, and AdaptiveKernel (the default). See the documentation for more information about these choices.

source
DelaunayTriangulation.has_acute_neighbouring_anglesMethod
has_acute_neighbouring_angles(predicates::AbstractPredicateKernel, enricher::BoundaryEnricher, i, j) -> Int, Vertex

Given a boundary edge (i, j), tests if the neighbouring angles are acute. The first returned value is the number of angles adjoining (i, j) that are acute (0, 1, or 2). The second returned value is the vertex that adjoins the edge (i, j) that is acute. If there are no such angles, or if there are two, then this returned vertex is 0.

(The purpose of this function is similar to segment_vertices_adjoin_other_segments_at_acute_angle.)

source
DelaunayTriangulation.split_subcurve!Function
split_subcurve!(enricher::BoundaryEnricher, i, j, predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Bool

Splits the curve associated with the edge (i, j) into two subcurves by inserting a point r between (i, j) such that the total variation of the subcurve is equal on (i, r) and (r, j). The returned value is a flag that is true if there was a precision issue, and false otherwise.

The predicate argument determines how predicates are computed, and should be one of ExactKernel, FastKernel, and AdaptiveKernel (the default). See the documentation for more information about these choices.

source
DelaunayTriangulation.test_visibilityMethod
test_visibility(predicates::AbstractPredicateKernel, enricher::BoundaryEnricher, i, j, k) -> Certificate

Tests if the vertex k is visible from the edge (i, j). Returns a Certificate which is

  • Invisible: If k is not visible from (i, j).
  • Visible: If k is visible from (i, j).

For this function, k should be inside the diametral circle of (i, j).

We say that k is invisibile from (i, j) if the edges (i, k) or (j, k) intersect any other boundary edges, or there is a hole between (i, j) and k.

Definition incompatibility

This is not the same definition used in defining constrained Delaunay triangulations, where visibility means visible from ANY point on the edge instead of only from the endpoints.

source
DelaunayTriangulation.triangulate_curve_boundedMethod
triangulate_curve_bounded(points::P;
segments=nothing,
boundary_nodes=nothing,
predicates::AbstractPredicateKernel=AdaptiveKernel(),
IntegerType::Type{I}=Int,
polygonise_n=4096,
coarse_n=0,
check_arguments=true,
delete_ghosts=false,
delete_empty_features=true,
recompute_representative_points=true,
rng::Random.AbstractRNG=Random.default_rng(),
insertion_order=nothing, 
kwargs...) where {P,I} -> Triangulation

Triangulates a curve-bounded domain defined by (points, segments, boundary_nodes). Please see triangulate for a description of the arguments. The only differences are:

  • insertion_order=nothing: This argument is ignored for curve-bounded domains.
  • polygonise_n=4096: For generating a high-resolution discretisation of a boundary initially for the construction of a PolygonHierarchy, many points are needed. This number of points is defined by polygonise_n, and must be a power of 2 (otherwise, the next highest power of 2 is used). See polygonise.
  • coarse_n=0: This is the number of points to use for initialising a curve-bounded domain via coarse_discretisation!. The default coarse_n=0 means the discretisation is performed until the maximum variation over any subcurve is less than π/2.
  • skip_points: This is still used, but it is ignored during the enrichment phase (see enrich_boundary!).

See also BoundaryEnricher and enrich_boundary!.

Refinement

To refine the mesh further beyond its initial coarse discretisation, as produced from this function, please see refine!.

source
DelaunayTriangulation.convert_boundary_curves!Method
convert_boundary_curves!(points, boundary_nodes, IntegerType) -> (NTuple{N, AbstractParametricCurve} where N, Vector)

Converts the provided points and boundary_nodes into a set of boundary curves and modified boundary nodes suitable for triangulation. In particular:

  1. The function gets boundary_curves from to_boundary_curves.
  2. boundary_nodes is replaced with a set of initial boundary nodes (from get_skeleton). These nodes come from evaluating each boundary curve at t = 0 and t = 1. In the case of a piecewise linear boundary, the vertices are copied directly. Note that not all control points of a CatmullRomSpline (which is_interpolating) will be added - only those at t = 0 and t = 1.
  3. The points are modified to include the new boundary nodes. If a point is already in points, it is not added again.

Arguments

  • points: The point set. This is modified in place with the new boundary points.
  • boundary_nodes: The boundary nodes to be converted. This is not modified in place.
  • IntegerType: The type of integer to use for the boundary nodes.

Output

  • boundary_curves: The boundary curves associated with boundary_nodes.
  • boundary_nodes: The modified boundary nodes.
source
DelaunayTriangulation.is_curve_boundedFunction
is_curve_bounded(tri::Triangulation) -> Bool 
is_curve_bounded(boundary_nodes) -> Bool

Returns true if tri is curve bounded, and false otherwise; similarly for the boundary_nodes method.

source
DelaunayTriangulation.test_visibilityFunction
test_visibility([kernel::AbstractPredicateKernel = AdaptiveKernel(),] tri::Triangulation, q, i) -> Certificate

Tests if the vertex i and the point q can see each other. Here, visibility means that the line segment joining the two does not intersect any segments.

Arguments

  • kernel::AbstractPredicateKernel = AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel (the default), and AdaptiveKernel. See the documentation for a further discussion of these methods.
  • tri: The Triangulation.
  • q: The point from which we are testing visibility.
  • i: The vertex we are testing visibility of.

Outputs

  • cert: A Certificate. This will be Visible if i is visible from q, and Invisible otherwise.
source
test_visibility([kernel::AbstractPredicateKernel=AdaptiveKernel(),] tri::Triangulation, u, v, i; shift=0.0, attractor=get_point(tri,i)) -> Certificate

Tests if the edge (u, v) and the point i can see each other. Here, visibility means that any point in the interior of (u, v) can see i. To test this, we only check 10 points equally spaced between u and v, excluding u and v.

Arguments

  • kernel::AbstractPredicateKernel=AdaptiveKernel(): How predicates are computed. See the documentation for information on the choices between FastKernel, ExactKernel, and AdaptiveKernel.
  • tri: The Triangulation.
  • u: The first vertex of the edge.
  • v: The second vertex of the edge.
  • i: The vertex we are testing visibility of.

Keyword Arguments

  • shift=0.0: The amount by which to shift each point on the edge towards attractor, i.e. if p is a point on the edge, then p .+ shift .* (attractor - p) is the point used to test visibility rather than p itself.
  • attractor=get_point(tri,i): Related to shift; see above.
  • predicates::AbstractPredicateKernel=AdaptiveKernel(): Method to use for computing predicates. Can be one of FastKernel, ExactKernel, and AdaptiveKernel. See the documentation for a further discussion of these methods.

Outputs

  • cert: A Certificate. This will be Visible if i is visible from (u, v), and Invisible otherwise.
source
test_visibility(predicates::AbstractPredicateKernel, enricher::BoundaryEnricher, i, j, k) -> Certificate

Tests if the vertex k is visible from the edge (i, j). Returns a Certificate which is

  • Invisible: If k is not visible from (i, j).
  • Visible: If k is visible from (i, j).

For this function, k should be inside the diametral circle of (i, j).

We say that k is invisibile from (i, j) if the edges (i, k) or (j, k) intersect any other boundary edges, or there is a hole between (i, j) and k.

Definition incompatibility

This is not the same definition used in defining constrained Delaunay triangulations, where visibility means visible from ANY point on the edge instead of only from the endpoints.

source