Data Structures

Here we list the internal data structures used throughout the package, along with functions related to working with the respective data structures.

MaxPriorityQueue

The MaxPriorityQueue is our implementation of a maximum priority queue; see, for example, the wikipedia page. The MaxPriorityQueue works similarly to the priority queue in DataStructures.jl. This struct is used for determining which triangles and edges to prioritise during mesh refinement, as well as which cells to prioritise during the quadtree construction for the pole of inaccessibility.

DelaunayTriangulation.MaxPriorityQueueType
MaxPriorityQueue{K, V}

Struct for a max priority queue.

Fields

  • data::Vector{Pair{K, V}}: The data of the queue, stored in a vector of key-value pairs mapping elements to their priority.
  • map::Dict{K, Int}: A dictionary mapping elements to their index in the data vector.
source
Base.firstMethod
first(queue::MaxPriorityQueue) -> Pair{K, V}

Returns the element with the highest priority in a queue, without removing it from queue.

source
Base.getindexMethod
getindex(queue::MaxPriorityQueue, key) 
queue[key]

Returns the priority of the element with key key in a queue.

source
Base.haskeyMethod
haskey(queue::MaxPriorityQueue, key) -> Bool

Returns true if the queue has an element with key key.

source
Base.isemptyMethod
isempty(queue::MaxPriorityQueue) -> Bool

Returns true if the queue is empty.

source
Base.lengthMethod
Base.length(queue::MaxPriorityQueue) -> Int

Returns the number of elements in queue.

source
Base.popfirst!Method
popfirst!(queue::MaxPriorityQueue{K, V}) where {K, V} -> Pair{K, V}

Removes and returns the element with the highest priority from the queue.

source
Base.push!Method
push!(queue::MaxPriorityQueue, pair)

Adds the key-value pair pair to the queue.

source
Base.setindex!Method
setindex!(queue::MaxPriorityQueue, priority, key)
queue[key] = priority

Sets the priority of the element with key key in a queue to priority, or adds the element to the queue if it is not already present.

source

Queue

The Queue is our implementation of a basic FIFO queue; see, for example, the wikipedia page. The implementation of a Queue in this package is very simple - it is just a vector. A block-based approach could also be implemented in the future, it just hasn't. The idea here is to avoid a dependency on DataStructures.jl. The Queue is used in a variety of places, such as during boundary enrichment in enrich_boundary! to determine which edges to process next.

DelaunayTriangulation.QueueType
Queue{T}

Struct for a first-in first-out queue.

Performance

Under the hood, Queue simply uses a Vector. This may not be as optimised compared to other implementations, e.g. DataStructure.jl's block-based approach with a Dequeue.

source
Base.eltypeMethod
eltype(queue::Queue{T}) -> Type{T}

Returns the type of elements stored in q.

source
Base.isemptyMethod
isempty(queue::Queue) -> Bool

Returns true if the queue is empty, false otherwise.

source
Base.lengthMethod
length(queue::Queue) -> Int

Returns the number of elements in the queue.

source
Base.popfirst!Method
popfirst!(queue::Queue)

Removes the element from the front of the queue and returns it.

source
Base.push!Method
push!(queue::Queue, item)

Adds item to the end of the queue.

source

BalancedBST

The BalancedBST is our implementation of a balanced binary search tree; see, for example, these notes. The BalancedBST is not currently used anywhere. Originally, it was intended for use with the Bentley-Ottman algorithm, but the complete algorithm is yet to be implemented; this algorithm was to be implemented to allow for users to provide intersecting segments that, through the algorithm, could be automatically split for the user. The BalancedBST implementation is fully functional for this purpose, the work in getting the actual algorithm up and going just has to be done.

DelaunayTriangulation.BalancedBSTType
mutable struct BalancedBST{K}

Struct representing a balanced binary search tree.

Fields

  • root::Union{Nothing,BalancedBSTNode{K}}: The root of the tree.
  • count::Int32: The number of nodes in the tree.
Duplicate keys

Nodes with duplicate keys are not supported. If a duplicate key is inserted, the tree will not be modified.

source
DelaunayTriangulation.BalancedBSTNodeType
mutable struct BalancedBSTNode{K}

Struct representing a node in a balanced binary search tree.

Fields

  • key::K: The key associated with the node.
  • height::Int8: The height of the node.
  • count::Int32: The number of nodes in the subtree rooted at this node, including this node.
  • parent::Union{Nothing,BalancedBSTNode{K}}: The parent of the node.
  • left::Union{Nothing,BalancedBSTNode{K}}: The left child of the node.
  • right::Union{Nothing,BalancedBSTNode{K}}: The right child of the node.

Constructor

BalancedBSTNode(key::K) where {K}

Constructs a new node with key key, height 1, and count 1. The parent, left child, and right child are set to nothing.

source
Base.delete!Method
delete!(tree::BalancedBST{K}, key::K) -> BalancedBST{K}

Deletes the node in tree with key key if it exists. Returns tree.

source
Base.findfirstMethod
findfirst(tree::BalancedBST{K}, key::K) -> Union{Nothing,BalancedBSTNode{K}}

Returns the node in tree with key key. If no such node exists, returns nothing.

source
Base.haskeyMethod
haskey(tree::BalancedBST{K}, key::K) -> Bool

Returns true if tree has a node with key key, false otherwise.

source
Base.minimumMethod
minimum(node::Union{BalancedBSTNode,Nothing}) -> Union{BalancedBSTNode,Nothing}

Returns the node with the minimum key in the subtree rooted at node. If node is nothing, returns nothing.

source
Base.push!Method
push!(tree::BalancedBST{K}, key::K)

Inserts key into tree if it is not already present.

source
DelaunayTriangulation.compute_balanceMethod
compute_balance(node::Union{Nothing,BalancedBSTNode{K}}) -> Int8

Computes the balance of the subtree rooted at node. This is the difference between the left and right heights.

source
DelaunayTriangulation.compute_countMethod
compute_count(node::Union{Nothing,BalancedBSTNode{K}}) -> Int32

Computes the count of the subtree rooted at node, i.e. the number of nodes in the subtree rooted at node, including node.

source
DelaunayTriangulation.delete_node!Method
delete_node!(node::Union{Nothing,BalancedBSTNode{K}}, key::K) -> Union{Nothing,BalancedBSTNode{K}}

Deletes the node with key key from the subtree rooted at node if it exists. Returns the new root of the subtree.

source
DelaunayTriangulation.get_leftMethod
get_left(node::BalancedBSTNode) -> Union{Nothing, BalancedBSTNode}

Returns the left child of node. If the node is nothing, returns nothing.

source
DelaunayTriangulation.get_rightMethod
get_right(node::BalancedBSTNode) -> Union{Nothing, BalancedBSTNode}

Returns the right child of node. If the node is nothing, returns nothing.

source
DelaunayTriangulation.insert_node!Method
insert_node!(node::Union{Nothing,BalancedBSTNode{K}}, key::K) -> Union{Nothing,BalancedBSTNode{K}}

Inserts key into the subtree rooted at node if it is not already present. Returns the new root of the subtree.

source
DelaunayTriangulation.rotate_left!Method
rotate_left!(parent::BalancedBSTNode{K}) -> BalancedBSTNode{K}

Rotates a subtree rooted at parent to the left, returning the new root of the subtree. This local operation is used to preserve the binary search tree property after inserting or deleting a node.

source
DelaunayTriangulation.rotate_right!Method
rotate_right!(parent::BalancedBSTNode{K}) -> BalancedBSTNode{K}

Rotates a subtree rooted at parent to the right, returning the new root of the subtree. This local operation is used to preserve the binary search tree property after inserting or deleting a node.

source

RTree

The RTree is our implementation of a simple RTree data structure, following heavily the original code in SpatialIndexing.jl. This tree is needed for triangulating curve-bounded domains. In particular, since the enrichment phase of triangulating a curve-bounded domain has no triangulation data structure that can be used for finding intersections between an edge's diametral circle and another vertex, the RTree is used to find these intersections.

DelaunayTriangulation.RTreeType
mutable struct RTree

Type for representing an R-tree with linear splitting.

Fields

  • root::Union{Branch,Leaf{Branch}}: The root of the R-tree.
  • num_elements::Int: The number of elements in the R-tree.
  • branch_cache::BranchCache: The cache of branch nodes.
  • twig_cache::TwigCache: The cache of twig nodes.
  • leaf_cache::LeafCache: The cache of leaf nodes.
  • fill_factor::Float64: The fill factor of the R-tree, i.e. the percentage of a node's capacity that should be filled after splitting.
  • free_cache::BitVector: A bit vector for keeping track of which indices in detached_cache are free.
  • detached_cache::Vector{Union{Branch,Leaf{Branch}}}: A cache of detached nodes, i.e. nodes that have been split from the R-tree. This is used for deleting nodes.
  • intersection_cache::NTuple{2,IntersectionCache}: Cache used for identifying intersections. Each element of the Tuple is its own cache, allowing for up to two intersection queries to be performed simultaneously. Note that this makes the R-tree non-thread-safe, and even non-safe when considering three or more intersection queries simultaneously.

Constructor

    RTree(; size_limit=100, fill_factor=0.7)

The size_limit is the node capacity. All node types have the same capacity.

source
DelaunayTriangulation.BoundaryRTreeType
BoundaryRTree{P}

Type for representing an R-tree of a boundary with an associated point set. The rectangular elements of the R-tree are the bounding box of the diametral circles of the boundary edges.

Fields

  • tree::RTree: The R-tree.
  • points::P: The point set.

Constructors

BoundaryRTree(points)
source
DelaunayTriangulation.BoundingBoxType
BoundingBox <: AbstractBoundingShape

Type for representing an axis-aligned bounding box, represented as a pair of interval I and J so that the bounding box is I × J.

Fields

  • x::BoundingInterval: The interval for the x-axis.
  • y::BoundingInterval: The interval for the y-axis.

Constructors

BoundingBox(x::BoundingInterval, y::BoundingInterval)
BoundingBox(a::Float64, b::Float64, c::Float64, d::Float64) = BoundingBox(BoundingInterval(a, b), BoundingInterval(c, d))
BoundingBox(p::NTuple{2,<:Number}) = BoundingBox(p[1], p[1], p[2], p[2])
source
DelaunayTriangulation.BranchType
mutable struct Branch <: AbstractNode

Type for representing a branch node in an R-tree.

Fields

  • parent::Union{Branch, Nothing}: The parent of the branch node.
  • bounding_box::BoundingBox: The bounding box of the branch node.
  • children::Union{Vector{Branch},Vector{Leaf{Branch}}}: The children of the branch node.
  • level::Int: The level of the branch node.

Constructor

Branch(parent::Union{Branch,Nothing}=nothing, ::Type{C}=Branch) where {C<:AbstractNode} = new(parent, InvalidBoundingBox, C[], 1)
source
DelaunayTriangulation.DiametralBoundingBoxType
DiametralBoundingBox

Type for representing a bounding box generated from an edge's diametral circle.

Fields

  • bounding_box::BoundingBox: The bounding box.
  • edge::NTuple{2,Int}: The generator edge.
source
DelaunayTriangulation.LeafType
mutable struct Leaf <: AbstractNode

Type for representing a leaf node in an R-tree.

Type parametrisation

Technically, this type should be referred to by Leaf{Branch}. Due to a lack of support for mutually recursive types or forward declarations, we have a parametric type in this struct's definition since Branch is not yet defined. In particular, Leaf is not a concrete type, whereas Leaf{Branch} is.

Fields

  • parent::Union{Branch, Nothing}: The parent of the leaf node.
  • bounding_box::BoundingBox: The bounding box of the leaf node.
  • children::Vector{DiametralBoundingBox}: The children of the leaf node.

Constructor

Leaf(parent::Union{Branch,Nothing}=nothing) = Leaf{Branch}(parent, InvalidBoundingBox, DiametralBoundingBox[])
source
DelaunayTriangulation.NodeCacheType
NodeCache{Node,Child}

Type for representing a cache of nodes whose children are of type Child. This is used for caching nodes that are detached from the R-tree, e.g. when a node is split.

Fields

  • cache::Vector{Node}: The cache of nodes.
  • size_limit::Int: The maximum number of nodes that can be cached.

Constructor

NodeCache{Node,Child}(size_limit::Int) where {Node,Child} = new{Node,Child}(Node[], size_limit)
source
DelaunayTriangulation.RTreeIntersectionCacheType
RTreeIntersectionCache

Type for representing a cache used for identifying intersections in an R-tree.

Fields

  • node_indices::Vector{Int}: A cache of indices used for identifying intersections.
  • need_tests::BitVector: A BitVector cache for keeping track of which indices in node_indices need to be tested for intersections.
source
DelaunayTriangulation.RTreeIntersectionIteratorType
RTreeIntersectionIterator

Type for representing an iterator over the elements in an R-tree that intersect with a bounding box.

Fields

  • tree::RTree: The R-tree.
  • bounding_box::BoundingBox: The bounding box to test for intersections with.
  • cache::RTreeIntersectionCache: The cache used for identifying intersections.
source
Base.:∩Method
intersect(r1::BoundingBox, r2::BoundingBox) -> BoundingBox
r1::BoundingBox ∩ r2::BoundingBox -> BoundingBox

Returns the intersection of r1 and r2. If the intersection is empty, returns InvalidBoundingBox.

source
Base.:∩Method
intersect(I::BoundingInterval, J::BoundingInterval) -> BoundingInterval
I::BoundingInterval ∩ J::BoundingInterval -> BoundingInterval

Returns the intersection of I and J. If the intersection is empty, returns InvalidBoundingInterval.

source
Base.:∪Method
union(r1::BoundingBox, r2::BoundingBox) -> BoundingBox
r1::BoundingBox ∪ r2::BoundingBox -> BoundingBox

Returns the union of r1 and r2, i.e. the smallest bounding box that contains both r1 and r2.

source
Base.:∪Method
union(I::BoundingInterval, J::BoundingInterval) -> BoundingInterval
I::BoundingInterval ∪ J::BoundingInterval -> BoundingInterval

Returns the union of I and J, combining their bounds; i.e. the smallest interval that contains both I and J.

source
Base.delete!Method
delete!(tree::BoundaryRTree, i, j)

Deletes the bounding box of the diametral circle of the edge between i and j in tree.

source
Base.inMethod
in(r1::BoundingBox, r2::BoundingBox) -> Bool
r1::BoundingBox ∈ r2::BoundingBox -> Bool

Tests whether r1 is in r2.

source
Base.inMethod
in(I::BoundingInterval, J::BoundingInterval) -> Bool
I::BoundingInterval ∈ J::BoundingInterval -> Bool

Tests whether the interval I is in the interval J.

source
Base.inMethod
in(a::Float64, I::BoundingInterval) -> Bool
a::Float64 ∈ I::BoundingInterval -> Bool

Tests whether a is in I.

source
Base.inMethod
in(p::NTuple{2,<:Number}, r::BoundingBox) -> Bool
p::NTuple{2,<:Number} ∈ r::BoundingBox -> Bool

Tests whether p is in r.

source
Base.insert!Method
insert!(tree::BoundaryRTree, i, j) -> Bool

Inserts the diametral circle of the edge between i and j into tree. Returns true if the tree's bounding boxes had to be adjusted and false otherwise.

source
Base.isemptyMethod
isempty(r::BoundingBox) -> Bool

Returns true if r is empty, i.e. if r.x or r.y is empty.

source
Base.isemptyMethod
isempty(I::BoundingInterval) -> Bool

Returns true if I is empty, i.e. if I.a or I.b is NaN or if length(I) < 0.

source
Base.isemptyMethod
isempty(cache::NodeCache) -> Bool

Returns true if cache is empty.

source
Base.lengthMethod
length(I::BoundingInterval) -> Float64

Returns the length of the interval I.

source
Base.lengthMethod
length(cache::NodeCache) -> Int

Returns the number of nodes in cache.

source
Base.pop!Method
pop!(cache::NodeCache) -> Node

Removes and returns the last node in cache.

source
Base.setindex!Method
setindex!(branch::Branch, child, i::Integer)
branch[i] = child

Sets the ith child of branch to be child, also updating the parent of child to be branch.

source
DelaunayTriangulation.bounding_boxMethod
bounding_box(points, i, j) -> DiametralBoundingBox

Returns the bounding box of the diametral circle of the points points[i] and points[j] with generator edge (i, j), returned as an DiametralBoundingBox.

source
DelaunayTriangulation.get_intersectionsMethod
get_intersections(tree::BoundaryRTree, i, j, k; cache_id=1) -> RTreeIntersectionIterator

Returns an RTreeIntersectionIterator over the elements in tree that intersect with the bounding box of the triangle (i, j, k). cache_id must be 1 or 2, and determines what cache to use for the intersection query.

source
DelaunayTriangulation.get_intersectionsMethod
get_intersections(tree::BoundaryRTree, i, j; cache_id=1) -> RTreeIntersectionIterator

Returns an RTreeIntersectionIterator over the elements in tree that intersect with the diametral circle of the edge between i and j. cache_id must be 1 or 2, and determines what cache to use for the intersection query.

source
DelaunayTriangulation.is_fullMethod
is_full(node::AbstractNode, tree::RTree) -> Bool

Returns true if node is full, i.e. if num_children(node) ≥ get_size_limit(tree).

source
DelaunayTriangulation.is_touchingMethod
is_touching(r1::BoundingBox, r2::BoundingBox) -> Bool

Tests whether r1 and r2 are touching, i.e. if they share a common boundary. This only considers interior touching.

source
DelaunayTriangulation.spawn_branch!Method
spawn_branch!(tree::RTree, bounding_box::BoundingBox, level) -> Branch

Returns a new branch node with bounding box bounding_box and level level from tree.

source
DelaunayTriangulation.spawn_node!Method
spawn_node!(tree::RTree, ::Type{N}, [bounding_box::BoundingBox], level) where {N} -> N

Returns a new node of type N with bounding box bounding_box and level level from tree. If bounding_box is not provided, it is replaced with InvalidBoundingBox.

source
DelaunayTriangulation.split_edge!Method
split_edge!(tree::BoundaryRTree, i, j, r)

Splits the diametral bounding box associated with (i, j) into two new boxes associated with the diametral circles of (i, r) and (j, r).

source
DelaunayTriangulation.QueryResultModule
QueryResult

An enum type for representing the result of an intersection query, with instances:

  • Contains: The bounding box contains the element.
  • Intersects: The bounding box intersects the element.
  • Outside: The bounding box is outside the element.
source

PolygonHierarchy

The PolygonHierarchy is a data structure used for representing hierarchical information about a set of polygons, such as the piecewise linear approximation of a curve-bounded domain's boundary. This structure is implemented as a collection of directed trees, where the disjoint trees each represent the disjoint parts of a domain, and the root of each tree contains the boundary curves of all the boundaries that are inside this tree.

DelaunayTriangulation.PolygonHierarchyType
PolygonHierarchy{I}

Struct used to define a polygon hierarchy. The hierarchy is represented as a forest of PolygonTrees.

Overlapping polygons

The polygons must not intersect any other polygon's boundaries.

Fields

  • polygon_orientations::BitVector: A BitVector of length n where n is the number of polygons in the hierarchy. The ith entry is true if the ith polygon is positively oriented, and false otherwise.
  • bounding_boxes::Vector{BoundingBox}: A Vector of BoundingBoxs of length n where n is the number of polygons in the hierarchy. The ith entry is the BoundingBox of the ith polygon.
  • trees::Dict{I,PolygonTree{I}}: A Dict mapping the index of a polygon to its PolygonTree. The keys of trees are the roots of each individual tree, i.e. the outer-most polygons.
  • reorder_cache::Vector{PolygonTree{I}}: A Vector used for caching trees to be deleted in [reorder_subtree!`](@ref).
One-based indexing

Note that the vector definitions for poylgon_orientations and bounding_boxes are treating the curves with the assumption that they are enumerated in the order 1, 2, 3, ....

Constructor

PolygonHierarchy{I}() where {I}

Constructs a PolygonHierarchy with no polygons.

source
DelaunayTriangulation.PolygonTreeType
mutable struct PolygonTree{I}

A tree structure used to define a polygon hierarchy.

Fields

  • parent::Union{Nothing,PolygonTree{I}}: The parent of the tree. If nothing, then the tree is the root.
  • children::Set{PolygonTree{I}}: The children of the tree.
  • index::I: The index of the tree. This is the index associated with the polygon.
  • height::Int: The height of the tree. This is the number of polygons in the tree that index is inside of. The root has height 0.

Constructor

PolygonTree{I}(parent::Union{Nothing,PolygonTree{I}}, index, height) where {I}

Constructs a PolygonTree with parent, index, and height, and no children.

source
DelaunayTriangulation.construct_polygon_hierarchyMethod
construct_polygon_hierarchy(points, boundary_nodes, boundary_curves; IntegerType=Int, n=4096) -> PolygonHierarchy{IntegerType}

Returns a PolygonHierarchy defining the polygon hierarchy for a given set of boundary_nodes that define a curve-bounded domain from the curves in boundary_curves. Uses polygonise to fill in the boundary curves.

Arguments

Keyword Arguments

  • IntegerType=Int: The integer type to use for indexing the polygons.
  • n=4096: The number of points to use for filling in the boundary curves in polygonise.
source
DelaunayTriangulation.find_treeMethod
find_tree(hierarchy::PolygonHierarchy, points, boundary_nodes, p) -> Union{Nothing,PolygonTree}

Finds a tree in hierarchy containing p.

Arguments

  • hierarchy::PolygonHierarchy: The PolygonHierarchy to search.
  • points: The point set.
  • boundary_nodes: The boundary nodes.
  • p: The point to test the trees of hierarchy against.

Output

  • nothing if p is not inside any tree in hierarchy, and the PolygonTree containing p otherwise.
source
DelaunayTriangulation.find_treeMethod
find_tree(hierarchy::PolygonHierarchy, points, boundary_nodes, tree::PolygonTree, p) -> PolygonTree

Finds a tree in hierarchy containing p that is a child of tree, assuming p is inside tree.

Arguments

  • hierarchy::PolygonHierarchy: The PolygonHierarchy to search.
  • points: The point set.
  • boundary_nodes: The boundary nodes.
  • tree::PolygonTree: The PolygonTree to search, assuming p is inside tree.
  • p: The point to test the children of tree against.

Output

  • tree if p is inside tree but none of its children, and a child containing p otherwise.
source
DelaunayTriangulation.is_in_treeMethod
is_in_tree(hierarchy::PolygonHierarchy, points, boundary_nodes, tree::PolygonTree, p) -> Bool

Tests if the point p is inside tree.

Arguments

  • hierarchy::PolygonHierarchy: The PolygonHierarchy containing tree.
  • points: The point set.
  • boundary_nodes: The boundary nodes.
  • tree::PolygonTree: The PolygonTree to test p against.
  • p: The point to test.

Output

  • true if p is inside tree, and false otherwise.
source
DelaunayTriangulation.reorder_hierarchy!Method
reorder_hierarchy!(hierarchy::PolygonHierarchy, points, boundary_nodes, new_tree::PolygonTree)

Given a new_tree that is not contained inside any other polygon in hierarchy, adds it into hierarchy. The existing trees are checked to see if they are contained inside new_tree, and if so, they are added as children of new_tree and removed from hierarchy.

Arguments

  • hierarchy::PolygonHierarchy: The PolygonHierarchy to add new_tree to.
  • points: The point set.
  • boundary_nodes: The boundary nodes.
  • new_tree::PolygonTree: The PolygonTree to add to hierarchy.
source
DelaunayTriangulation.reorder_subtree!Method
reorder_subtree!(hierarchy::PolygonHierarchy, points, boundary_nodes, tree::PolygonTree, new_tree)

Given a new_tree contained inside tree, adds it into hierarchy. The children of tree are reordered if necessary, in case they are now contained inside new_tree.

Arguments

  • hierarchy::PolygonHierarchy: The PolygonHierarchy to add new_tree to.
  • points: The point set.
  • boundary_nodes: The boundary nodes.
  • tree::PolygonTree: The PolygonTree to add new_tree to.
  • new_tree::PolygonTree: The PolygonTree to add to hierarchy and tree.
source
DelaunayTriangulation.set_bounding_box!Method
set_bounding_box!(hierarchy::PolygonHierarchy, index, bounding_box)

Sets the bounding box of the indexth polygon in hierarchy to bounding_box. If index is greater than the length of the bounding boxes vector, the vector is resized.

source
DelaunayTriangulation.set_orientation!Method
set_orientation!(hierarchy::PolygonHierarchy, index, orientation)

Sets the polygon orientation of the indexth polygon in hierarchy to orientation. If index is greater than the length of the polygon orientations vector, the vector is resized.

source
DelaunayTriangulation.set_tree!Method
set_tree!(hierarchy::PolygonHierarchy, index, tree)

Sets the PolygonTree of the indexth polygon in hierarchy to tree, or adds it if it is not an existing key. The index must be associated with an exterior polygon.

source

Adjacent

The Adjacent data structure is used for mapping edges to vertices (for triangulations) or polygons (for tessellations). The Adjacent structure itself is in the public API, as is get_adjacent.

DelaunayTriangulation.AdjacentType
Adjacent{IntegerType, EdgeType}

Struct for storing adjacency relationships for a triangulation.

Fields

  • adjacent::Dict{EdgeType, IntegerType}

The map taking edges (u, v) to w such that (u, v, w) is a positively oriented triangle in the underlying triangulation.

Constructors

Adjacent{IntegerType, EdgeType}()
Adjacent(adjacent::Dict{EdgeType, IntegerType})
source
DelaunayTriangulation.add_adjacent!Method
add_adjacent!(adj::Adjacent, uv, w)
add_adjacent!(adj::Adjacent, u, v, w)

Adds the adjacency relationship (u, v, w) to adj.

Examples

julia> using DelaunayTriangulation

julia> adj = DelaunayTriangulation.Adjacent{Int64, NTuple{2, Int64}}();

julia> DelaunayTriangulation.add_adjacent!(adj, 1, 2, 3)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 1 entry:
  (1, 2) => 3

julia> DelaunayTriangulation.add_adjacent!(adj, (2, 3), 1)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 2 entries:
  (1, 2) => 3
  (2, 3) => 1

julia> DelaunayTriangulation.add_adjacent!(adj, 3, 1, 2)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 3 entries:
  (1, 2) => 3
  (3, 1) => 2
  (2, 3) => 1
source
DelaunayTriangulation.add_triangle!Method
add_triangle!(adj::Adjacent, u, v, w)
add_triangle!(adj::Adjacent, T)

Adds the adjacency relationships defined from the triangle T = (u, v, w) to adj.

Examples

julia> using DelaunayTriangulation

julia> adj = DelaunayTriangulation.Adjacent{Int32, NTuple{2, Int32}}();

julia> add_triangle!(adj, 1, 2, 3)
Adjacent{Int32, Tuple{Int32, Int32}}, with map:
Dict{Tuple{Int32, Int32}, Int32} with 3 entries:
  (1, 2) => 3
  (3, 1) => 2
  (2, 3) => 1

julia> add_triangle!(adj, 6, -1, 7)
Adjacent{Int32, Tuple{Int32, Int32}}, with map:
Dict{Tuple{Int32, Int32}, Int32} with 6 entries:
  (1, 2)  => 3
  (3, 1)  => 2
  (6, -1) => 7
  (-1, 7) => 6
  (2, 3)  => 1
  (7, 6)  => -1
source
DelaunayTriangulation.delete_adjacent!Method
delete_adjacent!(adj::Adjacent, uv)
delete_adjacent!(adj::Adjacent, u, v)

Deletes the edge uv from adj.

Examples

julia> using DelaunayTriangulation

julia> adj = DelaunayTriangulation.Adjacent(Dict((2, 7) => 6, (7, 6) => 2, (6, 2) => 2, (17, 3) => -1, (-1, 5) => 17, (5, 17) => -1));

julia> DelaunayTriangulation.delete_adjacent!(adj, 2, 7)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 5 entries:
  (-1, 5) => 17
  (17, 3) => -1
  (6, 2)  => 2
  (5, 17) => -1
  (7, 6)  => 2

julia> DelaunayTriangulation.delete_adjacent!(adj, (6, 2))
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 4 entries:
  (-1, 5) => 17
  (17, 3) => -1
  (5, 17) => -1
  (7, 6)  => 2

julia> DelaunayTriangulation.delete_adjacent!(adj, 5, 17)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 3 entries:
  (-1, 5) => 17
  (17, 3) => -1
  (7, 6)  => 2
source
DelaunayTriangulation.delete_triangle!Method
delete_triangle!(adj::Adjacent, u, v, w)
delete_triangle!(adj::Adjacent, T)

Deletes the adjacency relationships defined from the triangle T = (u, v, w) from adj.

Examples

julia> using DelaunayTriangulation

julia> adj = DelaunayTriangulation.Adjacent{Int32, NTuple{2, Int32}}();

julia> add_triangle!(adj, 1, 6, 7);

julia> add_triangle!(adj, 17, 3, 5);

julia> adj
Adjacent{Int32, Tuple{Int32, Int32}}, with map:
Dict{Tuple{Int32, Int32}, Int32} with 6 entries:
  (17, 3) => 5
  (1, 6)  => 7
  (6, 7)  => 1
  (7, 1)  => 6
  (5, 17) => 3
  (3, 5)  => 17

julia> delete_triangle!(adj, 3, 5, 17)
Adjacent{Int32, Tuple{Int32, Int32}}, with map:
Dict{Tuple{Int32, Int32}, Int32} with 3 entries:
  (1, 6) => 7
  (6, 7) => 1
  (7, 1) => 6

julia> delete_triangle!(adj, 7, 1, 6)
Adjacent{Int32, Tuple{Int32, Int32}}, with map:
Dict{Tuple{Int32, Int32}, Int32}()
source
DelaunayTriangulation.get_adjacentMethod
get_adjacent(adj::Adjacent) -> Dict

Returns the adjacent map of adj.

Examples

julia> using DelaunayTriangulation

julia> d = Dict((1, 2) => 3, (2, 3) => 1, (3, 1) => 2);

julia> adj = DelaunayTriangulation.Adjacent(d)
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 3 entries:
  (1, 2) => 3
  (3, 1) => 2
  (2, 3) => 1

julia> get_adjacent(adj)
Dict{Tuple{Int64, Int64}, Int64} with 3 entries:
  (1, 2) => 3
  (3, 1) => 2
  (2, 3) => 1

julia> get_adjacent(adj) == d
true
source
DelaunayTriangulation.get_adjacentMethod
get_adjacent(adj::Adjacent{I, E}, uv::E) -> Vertex
get_adjacent(adj::Adjacent{I, E}, 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.

Examples

julia> using DelaunayTriangulation

julia> adj = DelaunayTriangulation.Adjacent(Dict((1, 2) => 3, (2, 3) => 1, (3, 1) => 2, (4, 5) => -1))
Adjacent{Int64, Tuple{Int64, Int64}}, with map:
Dict{Tuple{Int64, Int64}, Int64} with 4 entries:
  (4, 5) => -1
  (1, 2) => 3
  (3, 1) => 2
  (2, 3) => 1

julia> get_adjacent(adj, 4, 5)
-1

julia> get_adjacent(adj, (3, 1))
2

julia> get_adjacent(adj, (1, 2))
3

julia> get_adjacent(adj, 17, 5)
0

julia> get_adjacent(adj, (1, 6))
0
source

Adjacent2Vertex

The Adjacent2Vertex data structure is used for mapping vertices to edges for triangulations. The Adjacent2Vertex structure itself is in the public API, as is get_adjacent2vertex.

DelaunayTriangulation.Adjacent2VertexType
Adjacent2Vertex{IntegerType, EdgesType}

Struct for connectivity information about edges relative to vertices for a triangulation.

Fields

  • adjacent2vertex::Dict{IntegerType, EdgesType}

The map taking w to the set of all (u, v) such that (u, v, w) is a positively oriented triangle in the underlying triangle.

Constructors

Adjacent2Vertex{IntegerType, EdgesType}()
Adjacent2Vertex(adj2v::Dict{IntegerType, EdgesType})
source
DelaunayTriangulation.add_adjacent2vertex!Method
add_adjacent2vertex!(adj2v::Adjacent2Vertex, w, uv)
add_adjacent2vertex!(adj2v::Adjacent2Vertex, w, u, v)

Adds the edge uv to the set of edges E such that (u, v, w) is a positively oriented triangle in the underlying triangulation for each (u, v) ∈ E.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex{Int64, Set{NTuple{2, Int64}}}()
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}}()

julia> DelaunayTriangulation.add_adjacent2vertex!(adj2v, 1, (2, 3))
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 1 entry:
  1 => Set([(2, 3)])

julia> DelaunayTriangulation.add_adjacent2vertex!(adj2v, 1, 5, 7)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 1 entry:
  1 => Set([(5, 7), (2, 3)])

julia> DelaunayTriangulation.add_adjacent2vertex!(adj2v, 17, (5, -1))
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  17 => Set([(5, -1)])
  1  => Set([(5, 7), (2, 3)])
source
DelaunayTriangulation.add_triangle!Method
add_triangle!(adj2v::Adjacent2Vertex, u, v, w)
add_triangle!(adj2v::Adjacent2Vertex, T)

Adds the relationships defined by the triangle T = (u, v, w) into adj2v.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex{Int32, Set{NTuple{2, Int32}}}()
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}}()

julia> add_triangle!(adj2v, 17, 5, 8)
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}} with 3 entries:
  5  => Set([(8, 17)])
  8  => Set([(17, 5)])
  17 => Set([(5, 8)])

julia> add_triangle!(adj2v, 1, 5, 13)
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}} with 5 entries:
  5  => Set([(8, 17), (13, 1)])
  13 => Set([(1, 5)])
  8  => Set([(17, 5)])
  17 => Set([(5, 8)])
  1  => Set([(5, 13)])
source
DelaunayTriangulation.clear_empty_keys!Method
clear_empty_keys!(adj2v::Adjacent2Vertex)

Deletes all vertices w from adj2v such that get_adjacent2vertex(adj2v, w) is empty.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex{Int64, Set{NTuple{2, Int64}}}()
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}}()

julia> add_triangle!(adj2v, 1, 2, 3)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 3 entries:
  2 => Set([(3, 1)])
  3 => Set([(1, 2)])
  1 => Set([(2, 3)])

julia> delete_triangle!(adj2v, 2, 3, 1)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 3 entries:
  2 => Set()
  3 => Set()
  1 => Set()

julia> DelaunayTriangulation.clear_empty_keys!(adj2v)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}}()
source
DelaunayTriangulation.delete_adjacent2vertex!Method
delete_adjacent2vertex!(adj2v::Adjacent2Vertex, w, uv)
delete_adjacent2vertex!(adj2v::Adjacent2Vertex, w, u, v)

Deletes the edge uv from the set of edges E such that (u, v, w) is a positively oriented triangle in the underlying triangulation for each (u, v) ∈ E.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex(Dict(1 => Set(((2, 3), (5, 7), (8, 9))), 5 => Set(((1, 2), (7, 9), (8, 3)))))
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  5 => Set([(1, 2), (8, 3), (7, 9)])
  1 => Set([(8, 9), (5, 7), (2, 3)])

julia> DelaunayTriangulation.delete_adjacent2vertex!(adj2v, 5, 8, 3)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  5 => Set([(1, 2), (7, 9)])
  1 => Set([(8, 9), (5, 7), (2, 3)])

julia> DelaunayTriangulation.delete_adjacent2vertex!(adj2v, 1, (2, 3))
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  5 => Set([(1, 2), (7, 9)])
  1 => Set([(8, 9), (5, 7)])
source
DelaunayTriangulation.delete_adjacent2vertex!Method
delete_adjacent2vertex!(adj2v::Adjacent2Vertex, w)

Deletes the vertex w from adj2v.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex(Dict(1 => Set(((2, 3), (5, 7))), 5 => Set(((-1, 2), (2, 3)))))
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  5 => Set([(-1, 2), (2, 3)])
  1 => Set([(5, 7), (2, 3)])

julia> DelaunayTriangulation.delete_adjacent2vertex!(adj2v, 1)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 1 entry:
  5 => Set([(-1, 2), (2, 3)])

julia> DelaunayTriangulation.delete_adjacent2vertex!(adj2v, 5)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}}()
source
DelaunayTriangulation.delete_triangle!Method
delete_triangle!(adj2v::Adjacent2Vertex, u, v, w)
delete_triangle!(adj2v::Adjacent2Vertex, T)

Deletes the relationships defined by the triangle T =(u, v, w) from adj2v.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex{Int32, Set{NTuple{2, Int32}}}()
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}}()

julia> add_triangle!(adj2v, 1, 2, 3)
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}} with 3 entries:
  2 => Set([(3, 1)])
  3 => Set([(1, 2)])
  1 => Set([(2, 3)])

julia> add_triangle!(adj2v, 17, 5, 2)
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}} with 5 entries:
  5  => Set([(2, 17)])
  2  => Set([(3, 1), (17, 5)])
  17 => Set([(5, 2)])
  3  => Set([(1, 2)])
  1  => Set([(2, 3)])

julia> delete_triangle!(adj2v, 5, 2, 17)
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}} with 5 entries:
  5  => Set()
  2  => Set([(3, 1)])
  17 => Set()
  3  => Set([(1, 2)])
  1  => Set([(2, 3)])

julia> delete_triangle!(adj2v, 2, 3, 1)
Adjacent2Vertex{Int32, Set{Tuple{Int32, Int32}}} with map:
Dict{Int32, Set{Tuple{Int32, Int32}}} with 5 entries:
  5  => Set()
  2  => Set()
  17 => Set()
  3  => Set()
  1  => Set()
source
DelaunayTriangulation.get_adjacent2vertexMethod
get_adjacent2vertex(adj2v::Adjacent2Vertex, w) -> Edges

Returns the set of edges E such that (u, v, w) is a positively oriented triangle in the underlying triangulation for each (u, v) ∈ E.

Examples

julia> using DelaunayTriangulation

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex(Dict(1 => Set(((2, 3), (5, 7), (8, 9))), 5 => Set(((1, 2), (7, 9), (8, 3)))))
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  5 => Set([(1, 2), (8, 3), (7, 9)])
  1 => Set([(8, 9), (5, 7), (2, 3)])

julia> get_adjacent2vertex(adj2v, 1)
Set{Tuple{Int64, Int64}} with 3 elements:
  (8, 9)
  (5, 7)
  (2, 3)

julia> get_adjacent2vertex(adj2v, 5)
Set{Tuple{Int64, Int64}} with 3 elements:
  (1, 2)
  (8, 3)
  (7, 9)
source
DelaunayTriangulation.get_adjacent2vertexMethod
get_adjacent2vertex(adj2v::Adjacent2Vertex) -> Dict

Returns the adjacent2vertex map of adj2v.

Examples

julia> using DelaunayTriangulation

julia> e1 = Set(((1, 2), (5, 3), (7, 8)));

julia> e2 = Set(((2, 3), (13, 5), (-1, 7)));

julia> d = Dict(9 => e1, 6 => e2);

julia> adj2v = DelaunayTriangulation.Adjacent2Vertex(d)
Adjacent2Vertex{Int64, Set{Tuple{Int64, Int64}}} with map:
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  6 => Set([(13, 5), (-1, 7), (2, 3)])
  9 => Set([(1, 2), (7, 8), (5, 3)])

julia> get_adjacent2vertex(adj2v)
Dict{Int64, Set{Tuple{Int64, Int64}}} with 2 entries:
  6 => Set([(13, 5), (-1, 7), (2, 3)])
  9 => Set([(1, 2), (7, 8), (5, 3)])

julia> get_adjacent2vertex(adj2v) == d
true
source

Graph

The Graph data structure is used for storing connectivity information about vertices in a triangulation. The Graph structure itself is in the public API, as is get_graph and get_neighbours.

DelaunayTriangulation.GraphType
Graph{IntegerType}

Struct for storing neighbourhood relationships between vertices in a triangulation. This is an undirected graph.

Fields

  • vertices::Set{IntegerType}

The set of vertices in the underlying triangulation.

  • edges::Set{NTuple{2, IntegerType}}

The set of edges in the underlying triangulation.

  • neighbours::Dict{IntegerType, Set{IntegerType}}

The map taking vertices u to the set of all v such that (u, v) is an edge in the underlying triangulation.

Constructors

Graph{IntegerType}()
Graph(vertices::Set{IntegerType}, edges::Set{NTuple{2, IntegerType}}, neighbours::Dict{IntegerType, Set{IntegerType}})
source

Curves

There are many data structures used to define the curves we provide in this package, all subtyping the AbstractParametricCurve type. This type, and its subtypes, are all in the public API with the exception of PiecewiseLinear.

DelaunayTriangulation.AbstractParametricCurveType
abstract type AbstractParametricCurve <: Function end

Abstract type for representing a parametric curve parametrised over 0 ≤ t ≤ 1. The curves represented by this abstract type should not be self-intersecting, with the exception of allowing for closed curves.

The structs that subtype this abstract type must implement are:

Functions that are defined for all AbstractParametricCurve subtypes are:

Efficiently computing the total variation

The curves in this package evaluate the total variation not by evaluating the integral itself, but by taking care of the changes in orientation in the curve to efficiently compute it. This is done by using the orientation markers of the curves, obtained using orientation_markers, that stored in the field orientation_markers of these curves. The function marked_total_variation is then used to evaluate it. You may like to consider using these functions for any curve you wish to implement yourself, using e.g. the BezierCurve struct's implementation as a reference.

source
DelaunayTriangulation.LineSegmentType
LineSegment <: AbstractParametricCurve

Curve for representing a line segment, parametrised over 0 ≤ t ≤ 1. This curve can be using line_segment(t) and returns a tuple (x, y) of the coordinates of the point on the curve at t.

Fields

  • first::NTuple{2,Float64}: The first point of the line segment.
  • last::NTuple{2,Float64}: The last point of the line segment.
  • length::Float64: The length of the line segment.

Constructor

You can construct a LineSegment using

LineSegment(first, last)
source
DelaunayTriangulation.CircularArcType
CircularArc <: AbstractParametricCurve

Curve for representing a circular arc, parametrised over 0 ≤ t ≤ 1. This curve can be evaluated using circular_arc(t) and returns a tuple (x, y) of the coordinates of the point on the curve at t.

Fields

  • center::NTuple{2,Float64}: The center of the arc.
  • radius::Float64: The radius of the arc.
  • start_angle::Float64: The angle of the initial point of the arc, in radians.
  • sector_angle::Float64: The angle of the sector of the arc, in radians. This is given by end_angle - start_angle, where end_angle is the angle at last, and so might be negative for negatively oriented arcs.
  • first::NTuple{2,Float64}: The first point of the arc.
  • last::NTuple{2,Float64}: The last point of the arc.
  • pqr::NTuple{3, NTuple{2, Float64}}: Three points on the circle through the arc. This is needed for point_position_relative_to_curve.
Orientation

The angles start_angle and end_angle should be setup such that start_angle > end_angle implies a positively oriented arc, and start_angle < end_angle implies a negatively oriented arc. Moreover, they must be in [0°, 2π°).

Constructor

You can construct a CircularArc using

CircularArc(first, last, center; positive=true)

It is up to you to ensure that first and last are equidistant from center - the radius used will be the distance between center and first. The positive keyword argument is used to determine if the arc is positively oriented or negatively oriented.

source
DelaunayTriangulation.EllipticalArcType
EllipticalArc <: AbstractParametricCurve

Curve for representing an elliptical arc, parametrised over 0 ≤ t ≤ 1. This curve can be evaluated using elliptical_arc(t) and returns a tuple (x, y) of the coordinates of the point on the curve at t.

Fields

  • center::NTuple{2,Float64}: The center of the ellipse.
  • horz_radius::Float64: The horizontal radius of the ellipse.
  • vert_radius::Float64: The vertical radius of the ellipse.
  • rotation_scales::NTuple{2,Float64}: If θ is the angle of rotation of the ellipse, then this is (sin(θ), cos(θ)).
  • start_angle::Float64: The angle of the initial point of the arc measured from center, in radians. This angle is measured from the center prior to rotating the ellipse.
  • sector_angle::Float64: The angle of the sector of the arc, in radians. This is given by end_angle - start_angle, where end_angle is the angle at last, and so might be negative for negatively oriented arcs.
  • first::NTuple{2,Float64}: The first point of the arc.
  • last::NTuple{2,Float64}: The last point of the arc.

Constructor

You can construct an EllipticalArc using

EllipticalArc(first, last, center, major_radius, minor_radius, rotation; positive=true)

where rotation is the angle of rotation of the ellipse, in degrees. The positive keyword argument is used to determine if the arc is positively oriented or negatively oriented.

source
DelaunayTriangulation.BezierCurveType
BezierCurve <: AbstractParametricCurve

Curve for representing a Bezier curve, parametrised over 0 ≤ t ≤ 1. This curve can be evaluated using bezier_curve(t) and returns a tuple (x, y) of the coordinates of the point on the curve at t.

A good reference on Bezier curves is this.

See also BSpline and CatmullRomSpline.

Loops

This curve is only tested on loop-free curves (and closed curves that otherwise have no self-intersections). It is not guaranteed to work on curves with loops, especially for finding the nearest point on the curve to a given point.

Interpolation

Remember that Bezier curves are not interpolation curves. They only go through the first and last control points, but not the intermediate ones. If you want an interpolation curve, use CatmullRomSpline.

Fields

  • control_points::Vector{NTuple{2,Float64}}: The control points of the Bezier curve. The curve goes through the first and last control points, but not the intermediate ones.
  • cache::Vector{NTuple{2,Float64}}: A cache of the points on the curve. This is used to speed up evaluation of the curve using de Casteljau's algorithm.
  • lookup_table::Vector{NTuple{2,Float64}}: A lookup table for the Bezier curve, used for finding the point on the curve closest to a given point. The ith entry of the lookup table corresponds to the t-value i / (length(lookup_table) - 1).
  • orientation_markers::Vector{Float64}: The orientation markers of the curve. These are defined so that the orientation of the curve is monotone between any two consecutive markers. The first and last markers are always 0 and 1, respectively. See orientation_markers.
Concurrency

The cache is not thread-safe, and so you should not evaluate this curve in parallel.

Constructor

You can construct a BezierCurve using

BezierCurve(control_points::Vector{NTuple{2,Float64}}; lookup_steps=5000, kwargs...)

The keyword argument lookup_steps=100 controls how many time points in [0, 1] are used for the lookup table. The kwargs... are keyword arguments passed to orientation_markers.

source
DelaunayTriangulation.BSplineType
BSpline <: AbstractParametricCurve

Curve for representing a BSpline, parametrised over 0 ≤ t ≤ 1. This curve can be evaluated using b_spline(t) and returns a tuple (x, y) of the coordinates of the point on the curve at t.

See also BezierCurve and CatmullRomSpline.

Our implementation of a BSpline is based on https://github.com/thibauts/b-spline.

Loops

This curve is only tested on loop-free curves (and closed curves that otherwise have no self-intersections). It is not guaranteed to work on curves with loops, especially for finding the nearest point on the curve to a given point.

Interpolation

Remember that B-spline curves are not interpolation curves. They only go through the first and last control points, but not the intermediate ones. For an interpolating spline, see CatmullRomSpline.

Fields

  • control_points::Vector{NTuple{2,Float64}}: The control points of the BSpline. The curve goes through the first and last control points, but not the intermediate ones.
  • knots::Vector{Int}: The knots of the BSpline. You should not modify or set this field directly (in particular, do not expect any support for non-uniform B-splines).
  • cache::Vector{NTuple{2,Float64}}: A cache of the points on the curve. This is used to speed up evaluation of the curve using de Boor's algorithm.
  • lookup_table::Vector{NTuple{2,Float64}}: A lookup table for the B-spline curve, used for finding the point on the curve closest to a given point. The ith entry of the lookup table corresponds to the t-value i / (length(lookup_table) - 1).
  • orientation_markers::Vector{Float64}: The orientation markers of the curve. These are defined so that the orientation of the curve is monotone between any two consecutive markers. The first and last markers are always 0 and 1, respectively. See orientation_markers.

Constructor

You can construct a BSpline using

BSpline(control_points::Vector{NTuple{2,Float64}}; degree=3, lookup_steps=5000, kwargs...)

The keyword argument lookup_steps is used to build the lookup table for the curve. Note that the default degree=3 corresponds to a cubic B-spline curve. The kwargs... are keyword arguments passed to orientation_markers.

source
DelaunayTriangulation.CatmullRomSplineType
CatmullRomSpline <: AbstractParametricCurve

Curve for representing a Catmull-Rom spline, parametrised over 0 ≤ t ≤ 1. This curve can be evaluated using catmull_rom_spline(t) and returns a tuple (x, y) of the coordinates of the point on the curve at t.

For information on these splines, see e.g. this article and this article. Additionally, this article lists some nice properties of these splines.

Loops

This curve is only tested on loop-free curves (and closed curves that otherwise have no self-intersections). It is not guaranteed to work on curves with loops, especially for finding the nearest point on the curve to a given point.

Extension

Typically, Catmull-Rom splines are defined on segments of four control points, and drawn between the two interior control points. This creates an issue in that the first and last control points will not be joined to the spline. To overcome this, we extend the spline to the left and right during the evaluation of a spline, using the fields left and right defined below. The rules used for extending these points come from CatmullRom.jl, which extrapolates based on a Thiele-like cubic polynomial.

Fields

  • control_points::Vector{NTuple{2,Float64}}: The control points of the Catmull-Rom spline. The curve goes through each point.
  • knots::Vector{Float64}: The parameter values of the Catmull-Rom spline. The ith entry of this vector corresponds to the t-value associated with the ith control point. With an alpha parameter α, these values are given by knots[i+1] = knots[i] + dist(control_points[i], control_points[i+1])^α, where knots[1] = 0, and the vector is the normalised by dividing by knots[end].
  • lookup_table::Vector{NTuple{2,Float64}}: A lookup table for the Catmull-Rom spline, used for finding the point on the curve closest to a given point. The ith entry of the lookup table corresponds to the t-value i / (length(lookup_table) - 1).
  • alpha::Float64: The alpha parameter of the Catmull-Rom spline. This controls the type of the parametrisation, where alpha = 0 corresponds to uniform parametrisation, alpha = 1/2 corresponds to centripetal parametrisation, and alpha = 1 corresponds to chordal parametrisation. Must be in [0, 1]. For reasons similar to what we describe for tension below, we only support alpha = 1/2 for now. (If you do really want to change it, use the _alpha keyword argument in the constructor.)
  • tension::Float64: The tension parameter of the Catmull-Rom spline. This controls the tightness of the spline, with tension = 0 being the least tight, and tension = 1 leading to straight lines between the control points. Must be in [0, 1]. You can not currently set this to anything except 0.0 due to numerical issues with boundary refinement. (For example, equivariation splits are not possible if tension=1 since the curve is piecewise linear in that case, and for tension very close to 1, the equivariation split is not always between the provided times. If you really want to change it, then you can use the _tension keyword argument in the constructor - but be warned that this may lead to numerical issues and potentially infinite loops.)
  • left::NTuple{2,Float64}: The left extension of the spline. This is used to evaluate the spline on the first segment.
  • right::NTuple{2,Float64}: The right extension of the spline. This is used to evaluate the spline on the last segment.
  • lengths::Vector{Float64}: The lengths of the individual segments of the spline.
  • segments::Vector{CatmullRomSplineSegment}: The individual segments of the spline.
  • orientation_markers::Vector{Float64}: The orientation markers of the curve. These are defined so that the orientation of the curve is monotone between any two consecutive markers. The first and last markers are always 0 and 1, respectively. See orientation_markers.

Constructor

To construct a CatmullRomSpline, use

CatmullRomSpline(control_points::Vector{NTuple{2,Float64}}; lookup_steps=5000, kwargs...)

The keyword argument lookup_steps is used to build the lookup table for the curve, with lookup_steps giving the number of time points in [0, 1] used for the lookup table. The kwargs... are keyword arguments passed to orientation_markers.

source
DelaunayTriangulation.CatmullRomSplineSegmentType
CatmullRomSplineSegment <: AbstractParametricCurve

A single segment of a Camtull-Rom spline, representing by a cubic polynomial. Note that evaluating this curve will only draw within the two interior control points of the spline.

Based on this article.

Fields

  • a::NTuple{2,Float64}: The coefficient on .
  • b::NTuple{2,Float64}: The coefficient on .
  • c::NTuple{2,Float64}: The coefficient on t.
  • d::NTuple{2,Float64}: The constant in the polynomial.
  • p₁::NTuple{2,Float64}: The second control point of the segment.
  • p₂::NTuple{2,Float64}: The third control point of the segment.

With these fields, the segment is parametrised over 0 ≤ t ≤ 1 by q(t), where

q(t) = at³ + bt² + ct + d,

and q(0) = p₁ and q(1) = p₂, where the segment is defined by four control points p₀, p₁, p₂, and p₃.

This struct is callable, returning the interpolated point (x, y) at t as a NTuple{2,Float64}.

Constructor

To construct this segment, use

catmull_rom_spline_segment(p₀, p₁, p₂, p₃, α, τ)

Here, p₀, p₁, p₂, and p₃ are the four points of the segment (not a, b, c, and d), and α and τ are the parameters of the spline. The parameter α controls the type of the parametrisation, where

  • α = 0: Uniform parametrisation.
  • α = 1/2: Centripetal parametrisation.
  • α = 1: Chordal parametrisation.

The parameter τ is the tension, and controls the tightness of the segment. τ = 0 is the least tight, while τ = 1 leads to straight lines between the control points. Both α and τ must be in [0, 1].

source
DelaunayTriangulation.PiecewiseLinearType
PiecewiseLinear <: AbstractParametricCurve

Struct for representing a piecewise linear curve. This curve should not be interacted with or constructed directly. It only exists so that it can be an AbstractParametricCurve. Instead, triangulations use this curve to know that its boundary_nodes field should be used instead.

Existing methods

This struct does have fields, namely points and boundary_nodes (and boundarynodes should be a contiguous section). These are only used so that we can use this struct in [`anglebetween](@ref) easily. In particular, we need to allow for evaluating this curve att=0and att=1, and similarly for differentiating the curve att=0and att=1. For this, we have defined, lettingLbe aPiecewiseLinearcurve,L(0)to return the first point on the curve, and the last point otherwise (meaningL(h)is constant forh > 0`), and similarly for differentiation. Do NOT rely on the implementation of these methods.

source
DelaunayTriangulation._get_interval_for_get_circle_intersectionMethod
_get_interval_for_get_circle_intersection(c::AbstractParametricCurve, t₁, t₂, r) -> (Float64, Float64, NTuple{2, Float64})

Given a circle centered at c(t₁) with radius r, finds an initial interval for get_circle_intersection to perform bisection on to find a point of intersection. The returned interval is (tᵢ, tⱼ), where tᵢ is the parameter value of the first point in the interval and tⱼ is the parameter value of the last point in the interval. (The interval does not have to be sorted.) The third returned value is p = c(t₁).

source
DelaunayTriangulation.angle_betweenMethod
angle_between(c₁::AbstractParametricCurve, c₂::AbstractParametricCurve) -> Float64

Given two curves c₁ and c₂ such that c₁(1) == c₂(0), returns the angle between the two curves, treating the interior of the curves as being left of both.

source
DelaunayTriangulation.angle_betweenMethod
angle_between(L₁::LineSegment, L₂::LineSegment) -> Float64

Returns the angle between L₁ and L₂, assuming that L₁.last == L₂.first (this is not checked). For consistency with If the segments are part of some domain, then the line segments should be oriented so that the interior is to the left of both segments.

source
DelaunayTriangulation.get_circle_intersectionMethod
get_circle_intersection(c::AbstractParametricCurve, t₁, t₂, r) -> (Float64, NTuple{2,Float64})

Given a circle centered at c(t₁) with radius r, finds the first intersection of the circle with the curve after t₁ and less than t₂. It is assumed that such an intersection exists. The returned value is (t, q), where t is the parameter value of the intersection and q is the point of intersection.

source
DelaunayTriangulation.get_closest_pointMethod
get_closest_point(b::AbstractParametricCurve p) -> (Float64, NTuple{2,Float64})

Returns the t-value and the associated point q on the curve b that is nearest to p using a binary search. The search is done until the binary search interval is smaller than 1e-12. This function will only work if the curve b has a lookup table.

Loops

This function is only tested on loop-free curves. It is not guaranteed to work on curves with loops. Moreover, for this function to be accurate, you want the lookup table in b to be sufficiently dense.

source
DelaunayTriangulation.get_equidistant_splitMethod
get_equidistant_split(c::AbstractParametricCurve, t₁, t₂) -> Float64

Returns a value of t such that the arc length along c from t₁ to t is equal to the arc length along c from t to t₂. Uses the bisection method to compute the t-value.

source
DelaunayTriangulation.get_equivariation_splitMethod
get_equivariation_split(c::AbstractParametricCurve, t₁, t₂) -> Float64, Float64

Returns a value of t such that the total variation of c from t₁ to t is equal to the total variation of c from t to t₂. Uses the bisection method to compute the t-value. Also returns the new total variation of the two pieces.

source
DelaunayTriangulation.horizontal_inflection_pointsMethod
horizontal_inflection_points(c::AbstractParametricCurve; steps=200, iters = 50, tol = 1e-5) -> Vector{Float64}

Returns points t such that x''(t) = 0 and 0 ≤ t ≤ 1, where x'' is the second derivative of the x-coordinate of c. This function uses Newton's method to find the roots of x''. Note that these are only technically inflection points if x'''(t) ≠ 0 at these points, but this is not checked.

High-degree curves

For curves of very high degree, such as Bezier curves with steps control points or greater, this function might fail to return all inflection points.

Arguments

  • c::AbstractParametricCurve: The curve to find the horizontal inflection points of.

Keyword Arguments

  • steps=200: The number of t-values to use for seeding Newton's method. In particular, Newton's method is run for each initial value in LinRange(0, 1, steps).
  • iters=50: The number of iterations to run Newton's method for.
  • tol=1e-5: The tolerance to use for uniquetol. Also used for deciding whether a root is a valid root, i.e. if abs(x''(t)) > tol for a found root t, then t is not a valid root and is rejected.

Output

  • t: All inflection points, given in sorted order.
source
DelaunayTriangulation.horizontal_turning_pointsMethod
horizontal_turning_points(c::AbstractParametricCurve; steps=200, iters = 50, tol = 1e-5) -> Vector{Float64}

Returns points t such that x'(t) = 0 and 0 ≤ t ≤ 1, where x' is the derivative of the x-coordinate of c. This function uses Newton's method to find the roots of x'.

High-degree curves

For curves of very high degree, such as Bezier curves with steps control points or greater, this function might fail to return all turning points.

Arguments

  • c::AbstractParametricCurve: The curve to find the horizontal turning points of.

Keyword Arguments

  • steps=200: The number of t-values to use for seeding Newton's method. In particular, Newton's method is run for each initial value in LinRange(0, 1, steps).
  • iters=50: The number of iterations to run Newton's method for.
  • tol=1e-5: The tolerance to use for uniquetol. Also used for deciding whether a root is a valid root, i.e. if abs(x'(t)) > tol for a found root t, then t is not a valid root and is rejected.

Output

  • t: All turning points, given in sorted order.
source
DelaunayTriangulation.inflection_pointsMethod
inflection_points(c::AbstractParametricCurve; steps=200, iters = 50, tol = 1e-5) -> Vector{Float64}

Returns points t such that κ(t) = 0 and 0 ≤ t ≤ 1, where κ is the curvature of c. This function uses Newton's method to find the roots of κ.

High-degree curves

For curves of very high degree, such as Bezier curves with steps control points or greater, this function might fail to return all inflection points.

Arguments

  • c::AbstractParametricCurve: The curve to find the inflection points of.

Keyword Arguments

  • steps=200: The number of t-values to use for seeding Newton's method. In particular, Newton's method is run for each initial value in LinRange(0, 1, steps).
  • iters=50: The number of iterations to run Newton's method for.
  • tol=1e-5: The tolerance to use for uniquetol. Also used for deciding whether a root is a valid root, i.e. if abs(κ(t)) > tol for a found root t, then t is not a valid root and is rejected.
source
DelaunayTriangulation.orientation_markersMethod
orientation_markers(c::AbstractParametricCurve; steps=200, iters=50, tol=1e-5) -> Vector{Float64}

Finds all orientation markers of the AbstractParametricCurve c. These are points t where any of the following conditions hold (not necessarily simultaneously), letting c(t) = (x(t), y(t)):

  • x'(t) = 0
  • y'(t) = 0
  • κ(t; x) = 0, where κ(t; x) is the curvature of the component function x(t)
  • κ(t; y) = 0, where κ(t; y) is the curvature of the component function y(t)
  • κ(t) = 0, where κ is the curvature of c(t)

Note that the third and fourth conditions give all the inflection points of the component functions, and similarly for the fifth condition.

See also horizontal_turning_points, vertical_turning_points, horizontal_inflection_points, vertical_inflection_points, and inflection_points.

High-degree curves

For curves of very high degree, such as Bezier curves with steps control points or greater, this function might fail to return all inflection points.

Arguments

Keyword Arguments

  • steps=200: The number of equally spaced points to use for initialising Newton's method.
  • iters=50: How many iterations to use for Newton's method.
  • tol=1e-5: The tolerance used for determining if two t-values are the same.

Output

  • markers::Vector{Float64}: The t-values of the orientation markers of b. The returned vector is sorted, and also includes the endpoints 0 and 1; any t-values outside of [0, 1] are discarded, and multiplicity of any t is not considered (so the t-values in the returned vector are unique). These values can be used to split the curve into monotone pieces, meaning the orientation is monotone. These markers also guarantee that, over any monotone piece, the orientation changes by an angle of at most π/2.
source
DelaunayTriangulation.point_position_relative_to_curveMethod
point_position_relative_to_curve(e::AbstractParametricCurve, p) -> Certificate

Returns the position of the point p relative to the curve c. This function returns a [Certificate]:

  • Left: p is to the left of c.
  • Right: p is to the right of c.
  • On: p is on c.
source
DelaunayTriangulation.process_roots_and_residuals!Method
process_roots_and_residuals!(roots, residuals, tol) -> Vector{Float64}

Processes the roots and residuals of a root-finding algorithm. This function removes all NaN values from roots and residuals, sorts the roots in ascending order, and removes all roots with residuals greater than tol. The returned vector is the vector of roots with duplicates (i.e. roots that are within tol of each other) removed.

source
DelaunayTriangulation.protect_against_bad_division!Method
protect_against_bad_division!(roots, residuals, val, i) -> Bool

Protects against bad division in root-finding algorithms. This function checks if val is close to 0 or if roots[i] is outside of [0, 1]. If either of these conditions are true, then roots[i] and residuals[i] are set to NaN and true is returned. Otherwise, false is returned.

source
DelaunayTriangulation.vertical_inflection_pointsMethod
vertical_inflection_points(c::AbstractParametricCurve; steps=200, iters = 50, tol = 1e-5) -> Vector{Float64}

Returns points t such that y''(t) = 0 and 0 ≤ t ≤ 1, where y'' is the second derivative of the y-coordinate of c. This function uses Newton's method to find the roots of y''. Note that these are only technically inflection points if y'''(t) ≠ 0 at these points, but this is not checked.

High-degree curves

For curves of very high degree, such as Bezier curves with steps control points or greater, this function might fail to return all inflection points.

Arguments

  • c::AbstractParametricCurve: The curve to find the vertical inflection points of.

Keyword Arguments

  • steps=200: The number of t-values to use for seeding Newton's method. In particular, Newton's method is run for each initial value in LinRange(0, 1, steps).
  • iters=50: The number of iterations to run Newton's method for.
  • tol=1e-5: The tolerance to use for uniquetol. Also used for deciding whether a root is a valid root, i.e. if abs(y''(t)) > tol for a found root t, then t is not a valid root and is rejected.

Output

  • t: All inflection points, given in sorted order.
source
DelaunayTriangulation.vertical_turning_pointsMethod
vertical_turning_points(c::AbstractParametricCurve; steps=200, iters = 50, tol = 1e-5) -> Vector{Float64}

Returns points t such that y'(t) = 0 and 0 ≤ t ≤ 1, where y' is the derivative of the y-coordinate of c. This function uses Newton's method to find the roots of y'.

High-degree curves

For curves of very high degree, such as Bezier curves with steps control points or greater, this function might fail to return all turning points.

Arguments

  • c::AbstractParametricCurve: The curve to find the vertical turning points of.

Keyword Arguments

  • steps=200: The number of t-values to use for seeding Newton's method. In particular, Newton's method is run for each initial value in LinRange(0, 1, steps).
  • iters=50: The number of iterations to run Newton's method for.
  • tol=1e-5: The tolerance to use for uniquetol. Also used for deciding whether a root is a valid root, i.e. if abs(y'(t)) > tol for a found root t, then t is not a valid root and is rejected.

Output

  • t: All turning points, given in sorted order.
source

RepresentativeCoordinates

We use a RepresentativeCoordinates struct for storing the representative coordinates used for interpreting ghost vertices.

DelaunayTriangulation.RepresentativeCoordinatesType
RepresentativeCoordinates{IntegerType, NumberType}

A mutable struct for representing the coordinates of a representative point of polygon or a set of points.

Fields

  • x::NumberType: The x-coordinate of the representative point.
  • y::NumberType: The y-coordinate of the representative point.
  • n::IntegerType: The number of points represented by the representative point.
source

Cell

We use a Cell struct for representing an individual square in a quadtree during the computation of the pole of inaccessibility.

DelaunayTriangulation.CellType
Cell{T}

A cell in a grid. The cell is a square with side length 2half_width. The cell is centered at (x, y). The cell is assumed to live in a polygon.

Fields

  • x::T

The x-coordinate of the center of the cell.

  • y::T

The y-coordinate of the center of the cell.

  • half_width::T

The half-width of the cell.

  • dist::T

The distance from the center of the cell to the polygon.

  • max_dist::T

The maximum distance from the center of the cell to the polygon. This is dist + half_width * sqrt(2).

Constructors

Cell(x::T, y::T, half_width::T, points, boundary_nodes)

Constructs a cell with center (x, y) and half-width half_width. The cell is assumed to live in the polygon defined by points and boundary_nodes.

source

CellQueue

We use a CellQueue struct for storing the Cells in a quadtree during the computation of the pole of inaccessibility.

DelaunayTriangulation.CellQueueType
CellQueue{T}

A struct representing the priority queue of Cells, used for sorting the cells in a grid according to their maximum distance.

Fields

  • queue::MaxPriorityQueue{Cell{T},T}: The priority queue of cells, sorting according to maximum distance.

Constructors

CellQueue{T}()

Constructs a new CellQueue with elements of type Cell{T}.

source
Base.isemptyMethod
isempty(queue::CellQueue) -> Bool

Returns true if the queue is empty, and false otherwise.

source

ConvexHull

We use a ConvexHull struct to represent a convex hull. This struct is in the public API.

DelaunayTriangulation.ConvexHullType
ConvexHull{PointsType, IntegerType}

Struct for representing a convex hull. See also convex_hull.

Fields

  • points::PointsType: The underlying point set.
  • vertices::Vector{IntegerType}: The vertices of the convex hull, in counter-clockwise order. Defined so that vertices[begin] == vertices[end].

Constructors

ConvexHull(points, vertices)
convex_hull(points; IntegerType=Int)
source
DelaunayTriangulation.get_verticesMethod
get_vertices(convex_hull::ConvexHull) -> Vector{Vertices}

Returns the vertices of convex_hull. These are given in counter-clockwise order, and are defined so that the first and last vertices and equal.

source

Triangulation

We use a Triangulation struct to represent a triangulation. This struct is in the public API.

DelaunayTriangulation.TriangulationType
Triangulation{P,T,BN,W,I,E,Es,BC,BCT,BEM,GVM,GVR,BPL,C,BE}

Struct representing a triangulation, as constructed by triangulate.

Fields

  • points::P

The point set of the triangulation. Please note that this may not necessarily correspond to each point in the triangulation, e.g. some points may have been deleted - see each_solid_vertex for an iterator over each vertex in the triangulation.

  • triangles::T

The triangles in the triangulation. Each triangle is oriented counter-clockwise. If your triangulation has ghost triangles, some of these triangles will contain ghost vertices (i.e., vertices negative indices). Solid triangles can be iterated over using each_solid_triangle.

  • boundary_nodes::BN

The boundary nodes of the triangulation, if the triangulation is constrained; the assumed form of these boundary nodes is outlined in the docs. If your triangulation is unconstrained, then boundary_nodes will be empty and the boundary should instead be inspected using the convex hull field, or alternatively you can see lock_convex_hull!.

  • interior_segments::Es

Constrained segments appearing in the triangulation. These will only be those segments appearing off of the boundary. If your triangulation is unconstrained, then segments will be empty.

  • all_segments::Es

This is similar to segments, except this includes both the interior segments and the boundary segments. If your triangulation is unconstrained, then all_segments will be empty.

  • weights::W

The weights of the triangulation. If you are not using a weighted triangulation, this will be given by ZeroWeight(). Otherwise, the weights must be such that get_weight(weights, i) is the weight for the ith vertex. The weights should be Float64.

  • adjacent::Adjacent{I,E}

The Adjacent map of the triangulation. This maps edges (u, v) to vertices w such that (u, v, w) is a positively oriented triangle in triangles (up to rotation).

  • adjacent2vertex::Adjacent2Vertex{I,Es}

The Adjacent2Vertex map of the triangulation. This maps vertices w to sets S such that (u, v, w) is a positively oriented triangle in triangles (up to rotation) for all (u, v) ∈ S.

  • graph::Graph{I}

The Graph of the triangulation, represented as an undirected graph that defines all the neighbourhood information for the triangulation.

  • boundary_curves::BC

Functions defining the boundary curves of the triangulation, incase you are triangulating a curve-bounded domain. By default, this will be an empty Tuple, indicating that the boundary is as specified in boundary_nodes - a piecewise linear curve. If you are triangulating a curve-bounded domain, then these will be the parametric curves (see AbstractParametricCurve) you provided as a Tuple, where the ith element of the Tuple is associated with the ghost vertex -i, i.e. the ith section as indicated by ghost_vertex_map. If the ith boundary was left was a sequence of edges, then the function will be a PiecewiseLinear().

  • boundary_edge_map::BEM

This is a Dict from construct_boundary_edge_map that maps boundary edges (u, v) to their corresponding position in boundary_nodes.

  • ghost_vertex_map::GVM

This is a Dict that maps ghost vertices to their corresponding section in boundary_nodes, constructed by construct_ghost_vertex_map.

  • ghost_vertex_ranges::GVR

This is a Dict that maps ghost vertices to a range of all other ghost vertices that appear on the curve corresponding to the given ghost vertex, constructed by construct_ghost_vertex_ranges.

  • convex_hull::ConvexHull{P,I}

The ConvexHull of the triangulation, which is the convex hull of the point set points.

  • representative_point_list::BPL

The Dict of points giving RepresentativeCoordinates for each boundary curve, or for the convex hull if boundary_nodes is empty. These representative points are used for interpreting ghost vertices.

  • polygon_hierarchy::PolygonHierarchy{I}

The PolygonHierarchy of the boundary, defining the hierarchy of the boundary curves, giving information about which curves are contained in which other curves.

  • boundary_enricher::BE

The BoundaryEnricher used for triangulating a curve-bounded domain. If the domain is not curve-bounded, this is nothing.

  • cache::C

A TriangulationCache used as a cache for add_segment! which requires a separate Triangulation structure for use. This will not contain any segments or boundary nodes. Also stores segments useful for lock_convex_hull! and unlock_convex_hull!.

source
DelaunayTriangulation.TriangulationMethod
Triangulation(points, triangles, boundary_nodes; kwargs...) -> Triangulation

Returns the Triangulation corresponding to the triangulation of points with triangles and boundary_nodes.

Arguments

  • points: The points that the triangulation is of.
  • triangles: The triangles of the triangulation. These should be given in counter-clockwise order, with vertices corresponding to points. These should not include any ghost triangles.
  • boundary_nodes: The boundary nodes of the triangulation. These should match the specification given in the documentation or in check_args.

Keyword Arguments

  • IntegerType=Int: The integer type to use for the triangulation. This is used for representing vertices.
  • EdgeType=isnothing(segments) ? NTuple{2,IntegerType} : (edge_type ∘ typeof)(segments): The edge type to use for the triangulation.
  • TriangleType=NTuple{3,IntegerType}: The triangle type to use for the triangulation.
  • EdgesType=isnothing(segments) ? Set{EdgeType} : typeof(segments): The type to use for storing the edges of the triangulation.
  • TrianglesType=Set{TriangleType}: The type to use for storing the triangles of the triangulation.
  • weights=ZeroWeight(): The weights associated with the triangulation.
  • delete_ghosts=false: Whether to delete the ghost triangles after the triangulation is computed. This is done using delete_ghost_triangles!.

Output

source
DelaunayTriangulation.get_adjacentMethod
get_adjacent(tri::Triangulation) -> Adjacent

Returns the adjacency map of the triangulation. This is a map from each edge (u, v) to a vertex w such that (u, v, w) is a positively oriented triangle in tri.

See also Adjacent.

source
DelaunayTriangulation.get_all_segmentsMethod
get_all_segments(tri::Triangulation) -> Edges

Return all segments of the triangulation. This includes interior segments and boundary segments. Segments are edges that are forced to be in the triangulation.

source
DelaunayTriangulation.get_boundary_curvesMethod
get_boundary_curves(tri::Triangulation) -> NTuple{N, Function}

Returns the functions defining the boundaries of tri. If !is_curve_bounded(tri), then this returns an empty Tuple. Otherwise, this returns a Tuple of functions, one for each section of the boundary, where the ith element of the Tuple corresponds to the ith section of the boundary, which corresponds to the ghost vertex -i. For curves that are defined by boundary nodes rather than by a function, the function is PiecewiseLinear. For the other functions, these are all defined by t -> NTuple{2, Number}, where t ∈ [0, 1] and the NTuple{2, Number} is the coordinate on the curve at that t.

source
DelaunayTriangulation.get_boundary_edge_mapMethod
get_boundary_edge_map(tri::Triangulation) -> Dict

Returns the boundary edge map of the triangulation tri. This is a Dict that maps a boundary edge (u, v) to its position in get_boundary_nodes(tri). In particular, 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 (u, v) resides on, and u = get_boundary_nodes(boundary_nodes, index) and v = get_boundary_nodes(boundary_nodes, index + 1).

See also construct_boundary_edge_map.

source
DelaunayTriangulation.get_boundary_nodesMethod
get_boundary_nodes(tri::Triangulation) -> BoundaryNodes

Return the boundary nodes of the triangulation. This is only for triangulations with a constrained boundary. If the triangulation has no constrained boundary, then the boundary is instead given by its convex hull and this function returns an empty vector. See get_convex_hull.

source
DelaunayTriangulation.get_convex_hullMethod
get_convex_hull(tri::Triangulation) -> ConvexHull

Returns the convex hull of the points in tri. This is given as a ConvexHull object, where the vertices are sorted counter-clockwise and defined so that the first and last vertices are equal.

source
DelaunayTriangulation.get_ghost_vertex_mapMethod
get_ghost_vertex_map(tri::Triangulation) -> Dict

Returns the ghost vertex map of the triangulation tri. This is a Dict that maps ghost vertices to their associated section in boundary_nodes. There are three cases; below, I is integer_type(tri):

  • has_multiple_curves(tri)

Returns dict::Dict{I, NTuple{2, I}}, mapping ghost vertices i to Tuples (m, n) so that get_boundary_nodes(tri, m, n) are the boundary nodes associated with i, i.e. the nth section of the mth curve is associated with the ghost vertex i.

  • has_multiple_sections(tri)

Returns dict::Dict{I, I}, mapping ghost vertices i to n so that get_boundary_nodes(tri, n) are the boundary nodes associated with i, i.e. the nth section of the boundary is associated with the ghost vertex i.

  • otherwise

Returns dict::Dict{I, A}, mapping the ghost vertex i to get_boundary_nodes(tri), where A = typeof(get_boundary_nodes(tri)).

See also construct_ghost_vertex_map.

source
DelaunayTriangulation.get_interior_segmentsMethod
get_interior_segments(tri::Triangulation) -> Edges

Return the interior segments of the triangulation. These are segments that are forced to be in the triangulation - they are not the same as edges.

source
DelaunayTriangulation.get_trianglesMethod
get_triangles(tri::Triangulation) -> Triangles

Return the triangles of the triangulation. These triangles are all given in counter-clockwise order, and may include ghost triangles.

source

TriangulationCache

The TriangulationCache is what we store in the cache field of a triangulation.

DelaunayTriangulation.TriangulationCacheType
TriangulationCache{T,M,I,S}

A cache to be used as a field in Triangulation.

Fields

  • triangulation::T: The cached triangulation. This will only refer to an unconstrained triangulation, meaning it cannot have any segments or boundary nodes. It will contain the weights. This is used for constrained triangulations.
  • triangulation_2::T: An extra cached triangulation. This is needed for retriangulating fans for constrained triangulations.
  • marked_vertices::M: Marked vertices cache for use in constrained triangulations.
  • interior_segments_on_hull::I: Interior segments in the triangulation that also appear on the convex hull of tri. This is needed for lock_convex_hull! in case the convex hull also contains interior segments.
  • surrounding_polygon::S: The polygon surrounding the triangulation. This is needed for delete_point!.
  • fan_triangles::F: Triangles in a fan. This is needed for sorting fans for constrained triangulations.
Caches of caches

The triangulation cache itself does not have a cache. Instead, it stores a TriangulationCache(nothing).

Aliasing

The points of the cache's triangulation will be aliased to the points of the parent triangulation.

source
Base.empty!Method
empty!(cache::TriangulationCache)

Empties the cache by emptying the triangulation stored in it.

source

BoundaryEnricher

We use a BoundaryEnricher struct for storing information using during the enrichment phase of the triangulation of a curve-bounded domain.

DelaunayTriangulation.BoundaryEnricherType
BoundaryEnricher{P,B,C,I,BM,S}

Struct used for performing boundary enrichment on a curve-bounded boundary.

See also enrich_boundary!.

Fields

  • points::P: The point set.
  • boundary_nodes::B: The boundary nodes.
  • segments::S: The segments.
  • boundary_curves::C: The boundary curves.
  • polygon_hierarchy::PolygonHierarchy{I}: The polygon hierarchy.
  • parent_map::Dict{NTuple{2,I},I}: A map from an edge, represented as a Tuple, to the index of the parent curve in boundary_curves.
  • curve_index_map::Dict{I,I}: A map from a curve index to the index of the curve in boundary_curves.
  • boundary_edge_map::B: A map from a boundary node to the index of the curve in boundary_curves that it belongs to. See construct_boundary_edge_map.
  • spatial_tree::BoundaryRTree{P}: The BoundaryRTree used for spatial indexing.
  • queue::Queue{I}: A queue used for processing vertices during enrichment.
  • small_angle_complexes::Dict{I,Vector{SmallAngleComplex{I}}}: A map from an apex vertex to a list of all curves that define a small angle complex associated with that apex vertex.

The first three fields should be those associated with convert_boundary_curves!.

Constructor

BoundaryEnricher(points, boundary_nodes; IntegerType=Int, n=4096, coarse_n=0)

This constructor will use convert_boundary_curves! to convert points and boundary_nodes into a set of boundary curves and modified boundary nodes suitable for enrichment. The boundary nodes field will no longer aliased with the input boundary_nodes, although points will be. The polygon hierarchy is computed using construct_polygon_hierarchy. The argument n is used in polygonise for filling out the boundary temporarily in order to construct the PolygonHierarchy. The argument coarse_n defines the initial coarse discretisation through coarse_discretisation!; the default n=0 means that the coarse discretisation will be performed until the maximum total variation of a subcurve is less than π/2.

source
DelaunayTriangulation.SmallAngleComplexType
SmallAngleComplex{I}

Struct for representing a small-angle complex.

Fields

  • apex::I: The apex vertex of the complex.
  • members::Vector{SmallAngleComplexMember{I}}: The members of the complex.

Extended help

A small-angle complex is a set of curves that form a contiguous set of small angles, i.e. the angle between each consecutive pair of curves is less than 60°. The apex of the complex is the vertex that is shared by all of the curves.

source
DelaunayTriangulation.SmallAngleComplexMemberType
SmallAngleComplexMember{I}

Struct for representing a member of a small-angle complex.

Fields

  • parent_curve::I: The index of the parent curve in the boundary curves assoicated with the member. If this is 0, then this is instead a member of a complex around an interior segment.
  • next_edge::I: The next vertex after the apex in the boundary nodes associated with the member.
source
Base.append!Method
append!(complex::SmallAngleComplex, new_complex::SmallAngleComplex)

Appends the members of new_complex onto the members of complex.

source
Base.push!Method
push!(complex::SmallAngleComplex, member::SmallAngleComplexMember)

Pushes member onto the members of complex.

source
DelaunayTriangulation.construct_segment_mapMethod
construct_segment_map(segments, points, IntegerType) -> Dict{IntegerType, Vector{IntegerType}}

Returns the segment map of segments. This is a map that maps a vertex to all vertices that share a segment with that vertex. Segments are stored twice. The vertices associated with a vertex are sorted counter-clockwise, using the points argument to define the coordinates.

source
DelaunayTriangulation.construct_tree!Method
construct_tree!(enricher::BoundaryEnricher)

Constructs the spatial tree for enricher, modifying the spatial tree field in-place. The parent map must be correctly configured in order for this to be valid.

source
DelaunayTriangulation.get_boundary_curveMethod
get_boundary_curve(boundary_enricher::BoundaryEnricher, curve_index) -> AbstractParametricCurve

Returns the curve_indexth curve from the boundary curves in boundary_enricher.

source
DelaunayTriangulation.get_boundary_edge_mapMethod
get_boundary_edge_map(boundary_enricher::BoundaryEnricher, i, j)

Returns the value from the key (i, j) in the boundary edge map of boundary_enricher. The returned value is a Tuple (position, index) so that boundary_nodes = get_boundary_nodes(get_boundary_nodes(boundary_enricher), 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.get_circle_intersectionMethod
get_circle_intersection(enricher::BoundaryEnricher, curve_index, t₁, t₂, r) -> (Float64, NTuple{2,Float64})

Finds the intersection of the curve_indexth curve with the circle centered at the curve evaluated at t₁ with radius r. The argument t₂ defines the end of the subcurve to consider. The returned tuple is (t, p) where t is the parameter value of the intersection and p is the point of intersection.

source
DelaunayTriangulation.get_equivariation_splitMethod
get_equivariation_split(enricher::BoundaryEnricher, curve_index, t₁, t₂) -> Float64, Float64

Returns the equivariation split of the curve_indexth curve between t₁ and t₂. Also returns the total variation of the two pieces.

source
DelaunayTriangulation.get_parentMethod
get_parent(boundary_enricher::BoundaryEnricher{P,B,C,I}, i::I, j::I) -> I

Returns the parent of the edge (i, j) in boundary_enricher. If the edge is not in the parent map, then 0 is returned.

source
DelaunayTriangulation.get_small_angle_complexesFunction
get_small_angle_complexes(points, boundary_nodes, boundary_curves, segments=nothing; IntegerType=Int) -> Dict{IntegerType,Vector{SmallAngleComplex{IntegerType}}}

Returns a map from an apex vertex to a list of all curves that define a small angle complex associated with that apex vertex.

source
DelaunayTriangulation.is_segmentMethod
is_segment(enricher::BoundaryEnricher, i, j) -> Bool

Returns true if (i, j) or (j, i) is an interior segment of enricher, and false otherwise.

source
DelaunayTriangulation.is_small_angle_complex_memberMethod
is_small_angle_complex_member(enricher::BoundaryEnricher, i, j) -> Bool, I, IntegerType, IntegerType

Returns true if the edge (i, j) is a member of a small angle complex in enricher, and false otherwise.

Outputs

  • flag: true if the edge is a member of a small angle complex, and false otherwise.
  • apex: If the edge is a member of a small angle complex, then apex is the apex of the complex. Otherwise, apex is 0.
  • complex_id: If the edge is a member of a small angle complex, then complex_id is the index of the complex in the list of complexes associated with apex. Otherwise, complex_id is 0.
  • member_id: If the edge is a member of a small angle complex, then member_id is the index of the member in the list of members associated with complex_id. Otherwise, member_id is 0.
source
DelaunayTriangulation.partition_membersMethod
partition_members(complexes::Vector{SmallAngleComplex{I}}, points) where {I} -> Vector{SmallAngleComplex{I}}

Partitions the members of each complex in complexes into a new set of complexes. The complexes in complexes are assumed to be sorted in a counter-clockwise order around the apex of each complex. The partitioning is done so that the original set of members are now correctly split into their own complexes, since the original complexes might not have formed a properly contiguous set of small angles.

source
DelaunayTriangulation.polygoniseMethod
polygonise(points, boundary_nodes, boundary_curves; n=4096)

Fills out a set of points for a curve-bounded domain for use with PolygonHierarchy.

Warning

If the boundary curves are complicated so that they take a lot of points in order to be accurately resolved, then you should increase n.

Arguments

  • points: The point set.
  • boundary_nodes: The boundary nodes.
  • boundary_curves: The boundary curves.

Keyword Arguments

  • n=4096: The number of points to use for filling in each boundary curves.

Output

  • new_points: The points defining the filled out boundaries.
  • new_boundary_nodes: The boundary nodes associated with new_points.
Aliasing

If the boundary is not curve bounded, then new_points and new_boundary_nodes remain aliased with the input points and boundary_nodes.

source
DelaunayTriangulation.reorient_edgeMethod
reorient_edge(enricher::BoundaryEnricher, i, j) -> NTuple{2,Integer}

Given an edge (i, j), reorients it so that it is correctly oriented with the boundary. If (i, j) is instead an interior segment rather than a boundary edge, then (i, j) is returned.

source
DelaunayTriangulation.replace_next_edge!Method
replace_next_edge!(enricher::BoundaryEnricher, apex, complex_id, member_id, next_edge)

Replaces the next edge of the member_idth member of the complex_idth complex associated with apex with next_edge.

source
DelaunayTriangulation.replace_next_edgeMethod
replace_next_edge(member::SmallAngleComplexMember{I}, next_edge) where {I} -> SmallAngleComplexMember{I}

Returns a new SmallAngleComplexMember with the same parent curve as member but with next_edge as the next edge.

source
DelaunayTriangulation.split_boundary_edge!Function
split_boundary_edge!(enricher::BoundaryEnricher, i, j, r, update_boundary_nodes = Val(true))

Updates the fields of enricher after splitting a boundary edge (i, j) at the rth vertex. The update_boundary_nodes argument can be used to avoid inserting an additional boundary node when boundary_nodes was already updated somewhere else (e.g., we need this for mesh refinement which already updates the boundary_nodes which is aliased with the same field in the enricher).

source
DelaunayTriangulation.split_edge!Function
split_edge!(enricher::BoundaryEnricher, i, j, r, update_boundary_nodes = Val(true), update_segments = Val(true), is_interior = is_segment(enricher, i, j))

Updates the fields of enricher after splitting an edge (i, j) at the rth vertex. The update_boundary_nodes argument can be used to avoid inserting an additional boundary node when boundary_nodes was already updated somewhere else (e.g., we need this for mesh refinement which already updates the boundary_nodes which is aliased with the same field in the enricher). The same point goes for update_segments which can be used to avoid inserting an additional segment when segments was already updated somewhere else. The is_interior argument can be used to specify whether the edge is an interior segment or a boundary edge.

See also split_boundary_edge! and split_interior_segment!.

source
DelaunayTriangulation.split_interior_segment!Function
split_interior_segment!(enricher::BoundaryEnricher, i, j, r, update_segments = Val(true))

Updates the fields of enricher after splitting an interior segment (i, j) at the rth vertex. The update_segments argument can be used to avoid inserting an additional segment when segments was already updated somewhere else (e.g., we need this for mesh refinement which already updates the interior_segments which is aliased with the segments field in the enricher).

source

AbstractEach(Vertex/Edge/Triangle) Iterators

We use a variety of iterators for iterating over the vertices, edges, and triangles of a triangulation.

DelaunayTriangulation.EachGhostEdgeType
EachGhostEdge{E,T}

An iterator over all ghost edges in a triangulation.

Fields

  • edges::E: The iterator over all edges in the triangulation.
  • tri::T: The triangulation.
source
DelaunayTriangulation.EachGhostTriangleType
EachGhostTriangle{V,T}

An iterator over all ghost triangles in a triangulation.

Fields

  • triangles::V: The iterator over all triangles in the triangulation.
  • tri::T: The triangulation.
source
DelaunayTriangulation.EachGhostVertexType
EachGhostVertex{V,T}

An iterator over all ghost vertices in a triangulation.

Fields

  • vertices::V: The iterator over all vertices in the triangulation.
  • tri::T: The triangulation.
source
DelaunayTriangulation.EachSolidEdgeType
EachSolidEdge{E,T}

An iterator over all solid edges in a triangulation.

Fields

  • edges::E: The iterator over all edges in the triangulation.
  • tri::T: The triangulation.
source
DelaunayTriangulation.EachSolidTriangleType
EachSolidTriangle{V,T}

An iterator over all solid triangles in a triangulation.

Fields

  • triangles::V: The iterator over all triangles in the triangulation.
  • tri::T: The triangulation.
source
DelaunayTriangulation.EachSolidVertexType
EachSolidVertex{V,T}

An iterator over all solid vertices in a triangulation.

Fields

  • vertices::V: The iterator over all vertices in the triangulation.
  • tri::T: The triangulation.
source
DelaunayTriangulation.each_pointMethod
each_point(tri::Triangulation) -> Points

Returns an iterator over all points in tri.

Missing vertices

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 present in this iterator. It is recommended that you instead consider each_vertex, each_solid_vertex, or each_ghost_vertex together with get_point to get the coordinates.

source
DelaunayTriangulation.each_point_indexMethod
each_point_index(tri::Triangulation) -> Indices

Returns an iterator over all point indices in tri.

Missing vertices

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 point indices will still be present in this iterator. It is recommended that you instead consider each_vertex, each_solid_vertex, or each_ghost_vertex.

source

PointLocationHistory

We provide a means for storing the history of triangles encountered during point location, using a PointLocationHistory struct. The main motivation for this struct is for constrained triangulations.

DelaunayTriangulation.PointLocationHistoryType
PointLocationHistory{T,E,I}

History from using jump_and_march.

Fields

  • triangles::Vector{T}: The visited triangles.
  • collinear_segments::Vector{E}: Segments collinear with the original line pq using to jump.
  • collinear_point_indices::Vector{I}: This field contains indices to segments in collinear_segments that refer to points that were on the original segment, but there is no valid segment for them. We use manually fix this after the fact. For example, we could add an edge (1, 14), when really we mean something like (7, 14) which isn't a valid edge.
  • left_vertices::Vector{I}: Vertices from the visited triangles to the left of pq.
  • right_verices::Vector{I}: Vertices from the visited triangles to the right of pq.
source

IndividualTriangleStatistics

We provide an IndividualTriangleStatistics struct for storing statistics about individual triangles in a triangulation. This struct is in the public API, as listed in the API.

DelaunayTriangulation.IndividualTriangleStatisticsType
IndividualTriangleStatistics{T}

Struct storing statistics of a single triangle.

Fields

  • area::T: The area of the triangle.
  • lengths::NTuple{3,T}: The lengths of the edges of the triangle, given in sorted order.
  • circumcenter::NTuple{2,T}: The circumcenter of the triangle.
  • circumradius::T: The circumradius of the triangle.
  • angles::NTuple{3, T}: The angles of the triangle, given in sorted order.
  • radius_edge_ratio::T: The ratio of the circumradius to the shortest edge length.
  • edge_midpoints::NTuple{3,NTuple{2,T}}: The midpoints of the edges of the triangle.
  • aspect_ratio::T: The ratio of the inradius to the circumradius.
  • inradius::T: The inradius of the triangle.
  • perimeter::T: The perimeter of the triangle.
  • centroid::NTuple{2,T}: The centroid of the triangle.
  • offcenter::NTuple{2,T}: The offcenter of the triangle with radius-edge ratio cutoff β=1. See this paper.
  • sink::NTuple{2,T}: The sink of the triangle relative to the parent triangulation. See this paper.

Constructors

The constructor is

IndividualTriangleStatistics(p, q, r, sink = (NaN, NaN))

where p, q, and r are the coordinates of the triangle given in counter-clockwise order. sink is the triangle's sink. This must be provided separately since it is only computed relative to a triangulation, and so requires vertices rather than coordinates; see triangle_sink.

Extended help

The relevant functions used for computing these statistics are

source
DelaunayTriangulation.distance_to_offcenterMethod
distance_to_offcenter(β, ℓ) -> Number

Given a triangle with shortest edge length , computes the distance from the edge to the offcenter of the triangle with radius-edge ratio cutoff β.

source
DelaunayTriangulation.make_shortest_edge_firstMethod
make_shortest_edge_first(p, q, r, idx) -> (NTuple{2, Number}, NTuple{2, Number}, NTuple{2, Number})

Given a triangle (p, q, r), rotate it (preserving orientation) so that the shortest edge is first. The argument idx gives the index of the shortest edge, where idx == 1 means (p, q), idx == 2 means (q, r), and idx == 3 means (r, p).

source
DelaunayTriangulation.select_shortest_edge_for_offcenterMethod
select_shortest_edge_for_offcenter(p, q, r, c, ℓ²) -> (NTuple{2, Number}, NTuple{2, Number}, NTuple{2, Number})

Given a triangle (p, q, r) with more than one edge attaining the shortest length, selects the appropriate shortest edge for triangle_offcenter.

Arguments

  • p: The first vertex of the triangle.
  • q: The second vertex of the triangle.
  • r: The third vertex of the triangle.
  • c: The circumcenter of the triangle.
  • ℓ²: The squared shortest edge length.

The arguments should be so that (p, q, r) is positively oriented and ℓ² = |p - q|² is the squared shortest edge length.

Outputs

  • u: The first vertex of the rotated triangle.
  • v: The second vertex of the rotated triangle.
  • w: The third vertex of the rotated triangle.

These outputs (u, v, w) are a permutation of (p, q, r) (maintaining positive orientation) such that |m - c₁| is maximised over all other shortest edges, where m = (u + v)/2. If there is no unique maximiser, then the output is the permutation that is lexicographically smallest (i.e., sorted by x and then by y).

source
DelaunayTriangulation.squared_triangle_areaMethod
squared_triangle_area(ℓ₁², ℓ₂², ℓ₃²) -> Number

Compute the squared area of a triangle given the squares of its edge lengths. Heron's formula is used, so that the squared area is

\[A^2 = \dfrac{1}{16}\left[4\ell_1^2\ell_2^2 - \left(\ell_1^2 + \ell_2^2 - \ell_3^2\right)^2\right]..\]

See also squared_triangle_area_v2.

source
DelaunayTriangulation.squared_triangle_area_v2Method
squared_triangle_area_v2(ℓ₁², ℓ₂², ℓ₃²) -> Number

Compute the squared area of a triangle given the squares of its edge lengths, given in sorted order so that ℓ₁² ≤ ℓ₂² ≤ ℓ₃². This is a more numerically stable version of squared_triangle_area using the formula from Kahan (2014):

\[A^2 = \dfrac{1}{16}\left\{\left[\ell_3 + \left(\ell_2 + \ell_1\right)\right]\left[\ell_1 - \left(\ell_3 - \ell_2\right)\right]\left[\ell_1 + \left(\ell_3 - \ell_2\right)\right]\left[\ell_3 + \left(\ell_2 - \ell_1\right)\right]\right\}.\]

source
DelaunayTriangulation.squared_triangle_lengths_and_smallest_indexMethod
squared_triangle_lengths_and_smallest_index(p, q, r) -> (Number, Number, Number, Integer)

Computes the squared lengths of the edges of the triangle with coordinates p, q, r. The squared lengths are returned in sorted order, and the index of the shortest edge is returned as well. Here, the index refers to which edge in the order (p, q), (q, r), (q, p).

source
DelaunayTriangulation.triangle_anglesMethod
triangle_angles(p, q, r) -> (Number, Number, Number)

Computes the angles of a triangle with vertices p, q, and r. The formula for, say, the angle at p is given by

\[\theta_1 = \arctan\left(\dfrac{2A}{\left(p - q\right)\cdot\left(p - r\right)}\right),\]

where A is the area of the triangle. The angles are returned in sorted order.

source
DelaunayTriangulation.triangle_aspect_ratioMethod
triangle_aspect_ratio(inradius::Number, circumradius::Number) -> Number

Computes the aspect ratio of a triangle with inradius inradius and circumradius circumradius. The aspect ratio is given by

\[\tau = \dfrac{r_i}{r},\]

where $r_i$ is the inradius and $r$ is the circumradius.

source
DelaunayTriangulation.triangle_circumcenterFunction
triangle_circumcenter(p, q, r, A=triangle_area(p, q, r)) -> (Number, Number)

Computes the circumcenter of the triangle with coordinates (p, q, r). The circumcenter is given by

\[c_x = r_x + \dfrac{d_{11}d_{22} - d_{12}d_{21}}{4A}, \quad c_y = r_y + \dfrac{e_{11}e_{22} - e_{12}e_{21}}{4A},\]

where $d_{11} = \|p - r\|_2^2$, $d_{12} = p_y - r_y$, $d_{21} = \|q - r\|_2^2$, $d_{22} = q_y - r_y$, $e_{11} = p_x - r_x$ $e_{12} = d_{11}$, $e_{21} = q_x - r_x$, and $e_{22} = d_{21}$.

source
DelaunayTriangulation.triangle_circumradiusMethod
triangle_circumradius(A, ℓmin², ℓmed², ℓmax²) -> Number

Computes the circumradius of a triangle with area A and squared edge lengths ℓmin² ≤ ℓmed² ≤ ℓmax². The circumradius is given by

\[r = \dfrac{\ell_{\min}\ell_{\text{med}}\ell_{\max}}{4A}.\]

source
DelaunayTriangulation.triangle_inradiusMethod
triangle_inradius(A, perimeter) -> Number

Computes the inradius of a triangle with area A and perimeter perimeter. The inradius is given by

\[r_i = \dfrac{2A}{P},\]

where $P$ is the perimeter.

source
DelaunayTriangulation.triangle_lengthsMethod
triangle_lengths(p, q, r) -> (Number, Number, Number)

Computes the lengths of the edges of the triangle with coordinates p, q, r. The lengths are returned in sorted order.

source
DelaunayTriangulation.triangle_offcenterFunction
triangle_offcenter(p, q, r, c₁=triangle_circumcenter(p, q, r), β=1.0) -> (Number, Number)

Computes the off-center of the triangle (p, q, r).

Arguments

  • p, q, r: The coordinates of the triangle, given in counter-clockwise order.
  • c₁=triangle_circumcenter(p, q, r): The circumcenter of the triangle.
  • β=1.0: The radius-edge ratio cutoff.

Output

  • cx: The x-coordinate of the off-center.
  • cy: The y-coordinate of the off-center.
Difference in definitions

In the original this paper, the off-center is defined to instead be the circumcenter if it the triangle pqc₁ has radius-edge ratio less than β. Here, we just let the off-center be the point c so that pqc has radius-edge ratio of exactly β.

source
DelaunayTriangulation.triangle_perimeterMethod
triangle_perimeter(ℓmin::Number, ℓmed::Number, ℓmax::Number) -> Number

Computes the perimeter of a triangle with edge lengths ℓmin ≤ ℓmed ≤ ℓmax. The perimeter is given by

\[P = \ell_{\min} + \ell_{\text{med}} + \ell_{\max}.\]

source
DelaunayTriangulation.triangle_radius_edge_ratioMethod
triangle_radius_edge_ratio(circumradius::Number, ℓmin::Number) -> Number

Computes the radius-edge ratio of a triangle with circumradius circumradius and minimum edge length ℓmin, given by

\[\rho = \dfrac{r}{\ell_{\min}},\]

where $r$ is the circumradius and $\ell_{\min}$ is the shortest edge length.

source
DelaunayTriangulation.triangle_sinkFunction
triangle_sink(p, q, r, tri::Triangulation) -> (Number, Number)

Computes the sink of each triangle in tri. See this paper for more information.

Extended help

Sinks were introduced in this paper. For a given triangle T, the sink of T is defined as follows:

  1. If c, the circumcenter of T, is in the interior of T, then the sink of T is T.
  2. If T is a boundary triangle, then the sink of T is T.
  3. If neither 1 or 2, then the sink is defined as the sink of the triangle V, where V is the triangle adjoining the edge of T which intersects the line mc, where m is the centroid of T.

In cases where the triangulation has holes, this definition can lead to loops. In such a case, we just pick one of the triangles in the loop as the sink triangle.

source

TriangulationStatistics

We also provide a function for computing statistics about a triangulation, using the TriangulationStatistics struct. This struct is in the public API, as listed in the API.

DelaunayTriangulation.TriangulationStatisticsType
TriangulationStatistics{T,V,I}

A struct containing statistics about a triangulation.

Fields

  • num_vertices::I: The number of vertices in the triangulation.
  • num_solid_vertices::I: The number of solid vertices in the triangulation.
  • num_ghost_vertices::I: The number of ghost vertices in the triangulation.
  • num_edges::I: The number of edges in the triangulation.
  • num_solid_edges::I: The number of solid edges in the triangulation.
  • num_ghost_edges::I: The number of ghost edges in the triangulation.
  • num_triangles::I: The number of triangles in the triangulation.
  • num_solid_triangles::I: The number of solid triangles in the triangulation.
  • num_ghost_triangles::I: The number of ghost triangles in the triangulation.
  • num_boundary_segments::I: The number of boundary segments in the triangulation.
  • num_interior_segments::I: The number of interior segments in the triangulation.
  • num_segments::I: The number of segments in the triangulation.
  • num_convex_hull_vertices::I: The number of vertices on the convex hull of the triangulation.
  • smallest_angle::V: The smallest angle in the triangulation.
  • largest_angle::V: The largest angle in the triangulation.
  • smallest_area::V: The smallest area of a triangle in the triangulation.
  • largest_area::V: The largest area of a triangle in the triangulation.
  • smallest_radius_edge_ratio::V: The smallest radius-edge ratio of a triangle in the triangulation.
  • largest_radius_edge_ratio::V: The largest radius-edge ratio of a triangle in the triangulation.
  • area::V: The total area of the triangulation.
  • individual_statistics::Dict{T,IndividualTriangleStatistics{V}}: A map from triangles in the triangulation to their individual statistics. See IndividualTriangleStatistics.

Constructors

To construct these statistics, use statistics, which you call as statistics(tri::Triangulation).

source

InsertionEventHistory

For mesh refinement we need a way to identify what happens to a triangulation after a point is added, in case we need to reverse the insertion. For this, we use InsertionEventHistory internally.

Base.empty!Method
empty!(events::InsertionEventHistory)

Empties events by emptying all of its fields.

source
DelaunayTriangulation.each_added_segmentMethod
each_deleted_triangle(events::InsertionEventHistory) -> Iterator

Returns an iterator over the triangles that were deleted from the triangulation during the insertion of a point.

source
DelaunayTriangulation.split_boundary_edge!Method
split_boundary_edge!(events::InsertionEventHistory, u, v, new_point)

Add the edge (u, v) to the deleted_boundary_segments of events and add the edges (u, new_point) and (new_point, v) to the added_boundary_segments of events.

source
DelaunayTriangulation.undo_boundary_segment_changes!Method
undo_boundary_segment_changes!(tri::Triangulation, events::InsertionEventHistory)

Undoes any changes to the boundary segments in tri that were made after an insertion event, as recorded in events, assuming the inserted vertex is num_points(tri).

source
DelaunayTriangulation.undo_insertion!Function
undo_insertion!(tri::Triangulation, events::InsertionEventHistory, pop=Val(true))

Undoes the insertion of the most recent vertex into tri, assuming that its insertion history has been recorded into events and the vertex is num_points(tri).

If you do not want to delete the latest vertex from the triangulation, set pop to Val(false).

source

RefinementConstraints

For mesh refinement, internally we store the user-provided constraints inside their own struct RefinementConstraints.

DelaunayTriangulation.RefinementConstraintsType
RefinementConstraints{F}

A struct for storing constraints for mesh refinement.

Fields

  • min_angle=0.0: The minimum angle of a triangle.
  • max_angle=180.0: The maximum angle of a triangle.
  • min_area=0.0: The minimum area of a triangle.
  • max_area=Inf: The maximum area of a triangle.
  • max_radius_edge_ratio=csd(min_angle) / 2: The maximum radius-edge ratio of a triangle. This is computed from min_angle - you cannot provide a value for it yourself.
  • max_points=typemax(Int): The maximum number of vertices allowed in the triangulation. Note that this refers to num_solid_vertices, not the amount returned by num_points.
  • seiditous_angle=20.0: The inter-segment angle used to define seditious edges in degrees. Should not be substantially smaller than 20.0° or any greater than 60.0°.
  • custom_constraint::F=(tri, triangle) -> false: A custom constraint function. This should take a Triangulation and a triangle as arguments, and return true if the triangle violates the constraints and false otherwise.
source

RefinementQueue

The mesh refinement algorithm we use requires a method for prioritising which segments and triangles need to be split. Using a MaxPriorityQueue, so that large segments and large triangles are prioritised, we use a RefinementQueue internally to determine this ordering.

DelaunayTriangulation.RefinementQueueType
RefinementQueue{T,E,F}

Struct defining a pair of priority queues for encroached segments and poor-quality triangles.

Fields

  • segments::MaxPriorityQueue{E,F}: A priority queue of encroached segments, where the priorities are the squared edge lengths.
  • triangles::MaxPriorityQueue{T,F}: A priority queue of poor-quality triangles, where the priorities are the radius-edge ratios.

Constructor

The default constructor is available, but we also provide

RefinementQueue(tri::Triangulation)

which will initialise this struct with empty queues with the appropriate types.

source
Base.getindexMethod
getindex(queue::RefinementQueue{T,E,F}, triangle::T) -> F
queue[triangle] -> F

Return the radius-edge ratio of triangle in queue.

source
Base.haskeyMethod
haskey(queue::RefinementQueue{T,E,F}, segment::E) -> Bool

Return true if queue has segment or its reverse, and false otherwise.

source
Base.haskeyMethod
haskey(queue::RefinementQueue{T,E,F}, triangle::T) -> Bool

Return true if queue has triangle or any of its counter-clockwise rotations, and false otherwise.

source
Base.isemptyMethod
isempty(queue::RefinementQueue) -> Bool

Return true if queue has no segments or triangles, false otherwise.

source
Base.setindex!Method
Bassetindex!(queue::RefinementQueue{T,E,F}, ρ::F, triangle::T) where {T,E,F}
queue[triangle] = ρ

Add a triangle to queue whose radius-edge ratio is ρ. If the triangle is already in the queue, its priority is updated to ρ.

source
Base.setindex!Method
setindex!(queue::RefinementQueue{T,E,F}, ℓ²::F, segment::E) where {T,E,F}
queue[segment] = ℓ²

Add a segment to queue whose squared length is ℓ². If the segment is already in the queue, its priority is updated to .

source

RefinementArguments

There are many arguments that need to be passed around during mesh refinement. To simplify this, we use a RefinementArguments struct internally.

DelaunayTriangulation.RefinementArgumentsType
RefinementArguments{Q,C,H,I,E,T,R}

A struct for storing arguments for mesh refinement.

Fields

  • queue::Q: The RefinementQueue.
  • constraints::C: The RefinementConstraints.
  • events::H: The InsertionEventHistory.
  • min_steiner_vertex::I: The minimum vertex of a Steiner point. All vertices greater than or equal to this can be considered as Steiner vertices.
  • segment_list::Set{E}: The set of segments in the triangulation before refinement.
  • segment_vertices::Set{I}: Set of vertices that are vertices of segments in segment_list.
  • midpoint_split_list::Set{I}: Set of vertices that are centre-splits of encroached edges.
  • offcenter_split_list::Set{I}: Set of vertices that are off-centre splits of encroached edges.
  • use_circumcenter::Bool: Whether to use circumcenters for Steiner points, or the more general approach of Erten and Üngör (2009).
  • use_lens::Bool: Whether to use diametral lens (true) or diametral circles (false) for the defining encroachment.
  • steiner_scale::T: The factor by which to scale the Steiner points closer to the triangle's shortest edge.
  • locked_convex_hull::Bool: Whether the convex hull of the triangulation had to be locked for refinement.
  • had_ghosts::Bool: Whether the triangulation initially had ghost triangles or not.
  • rng::R: The random number generator.
  • concavity_protection::Bool: Whether to use concavity protection or not for jump_and_march. Most likely not needed, but may help in pathological cases.

Constructors

In addition to the default constructor, we provide

RefinementArguments(tri::Triangulation; kwargs...)

for constructing this struct. This constructor will lock the convex hull and add ghost triangles to tri if needed (refine! will undo these changes once the refinement is finished))

source
DelaunayTriangulation.RefinementArgumentsMethod
RefinementArguments(tri::Triangulation; kwargs...) -> RefinementArguments

Initialises the RefinementArguments for the given Triangulation, tri. The kwargs... match those from refine!.

Mutation

If tri has no ghost triangles, it will be mutated so that it has them. Similarly, if the triangulation has no constrained boundary, then the convex hull will be locked so that it is treated as a constrained boundary. These changes will be undone in refine! once the refinement is finished.

source
DelaunayTriangulation.is_freeMethod
is_free(args::RefinementArguments, u) -> Bool

Returns true if u is a free vertex, and false otherwise. A free vertex is a Steiner vertex (meaning a non-input vertex) that is not part of a segment or subsegment.

source
DelaunayTriangulation.keep_iteratingMethod
keep_iterating(tri::Triangulation, args::RefinementArguments) -> Bool

Returns true if the refinement should continue, and false otherwise. The check is based on whether the RefinementQueue is empty or not, and whether the number of points in the triangulation is less than or equal to the maximum number of points allowed by the RefinementConstraints.

source
DelaunayTriangulation.keep_splittingMethod
keep_splitting(tri::Triangulation, args::RefinementArguments) -> Bool

Returns true if more encroached segments need to be split, and false otherwise. The check is based on whether the segment queue in the RefinementQueue of args is empty or not, and whether the number of points in the triangulation is less than or equal to the maximum number of points allowed by the RefinementConstraints in args.

source

VoronoiTessellation

We use a VoronoiTessellation struct to represent a Voronoi tessellation. This struct is in the public API.

DelaunayTriangulation.VoronoiTessellationType
VoronoiTessellation{Tr<:Triangulation,P,I,T,S,E}

Struct for representing a Voronoi tessellation.

See also voronoi.

Fields

  • triangulation::Tr: The underlying triangulation. The tessellation is dual to this triangulation, although if the underlying triangulation is constrained then this is no longer the case (but it is still used).
  • generators::Dict{I,P}: A Dict that maps vertices of generators to coordinates. These are simply the points present in the triangulation. A Dict is needed in case the triangulation is missing some points.
  • polygon_points::Vector{P}: The points defining the coordinates of the polygons. The points are not guaranteed to be unique if a circumcenter appears on the boundary or you are considering a clipped tessellation. (See also get_polygon_coordinates.)
  • polygons::Dict{I,Vector{I}}: A Dict mapping polygon indices (which is the same as a generator vertex) to the vertices of a polygon. The polygons are given in counter-clockwise order and the first and last vertices are equal.
  • circumcenter_to_triangle::Dict{I,T}: A Dict mapping a circumcenter index to the triangle that contains it. The triangles are sorted such that the minimum vertex is last.
  • triangle_to_circumcenter::Dict{T,I}: A Dict mapping a triangle to its circumcenter index. The triangles are sorted such that the minimum vertex is last.
  • unbounded_polygons::Set{I}: A Set of indices of the unbounded polygons.
  • cocircular_circumcenters::S: A Set of indices of circumcenters that come from triangles that are cocircular with another triangle's vertices, and adjoin said triangles.
  • adjacent::Adjacent{I,E}: The adjacent map. This maps an oriented edge to the polygon that it belongs to.
  • boundary_polygons::Set{I}: A Set of indices of the polygons that are on the boundary of the tessellation. Only relevant for clipped tessellations, otherwise see unbounded_polygons.
source
DelaunayTriangulation.add_adjacent!Method
add_adjacent!(vor::VoronoiTessellation, ij, k)
add_adjacent!(vor::VoronoiTessellation, i, j, k)

Adds the adjacency relationship (i, j) ⟹ k between the oriented edge (i, j) and polygon index k to the Voronoi tessellation vor.

source
DelaunayTriangulation.add_polygon!Method
add_polygon!(vor::VoronoiTessellation, B, i)

Adds, or replaces, the polygon associated with the index i with B. B should be a counter-clockwise sequence of vertices, with B[begin] == B[end].

source
DelaunayTriangulation.get_adjacentMethod
get_adjacent(vor::VoronoiTessellation, ij) -> Index 
get_adjacent(vor::VoronoiTessellation, i, j) -> Index

Gets the polygon index associated with the oriented edge ij in the Voronoi tessellation vor.

source
DelaunayTriangulation.get_adjacentMethod
get_adjacent(vorn::VoronoiTessellation) -> Adjacent{Index,Edge}

Gets the adjacency information of the Voronoi tessellation vorn as an Adjacent object. This object maps oriented edges to the polygons that they belong to.

source
DelaunayTriangulation.get_circumcenter_to_triangleMethod
get_circumcenter_to_triangle(vorn::VoronoiTessellation) -> Dict{Index,Triangle}

Gets the circumcenters of the Voronoi tessellation vorn as a Dict, mapping circumcenter indices to their corresponding triangles. The triangles are sorted so that the minimum vertex is last.

source
DelaunayTriangulation.get_cocircular_circumcentersMethod
get_cocircular_circumcenters(vorn::VoronoiTessellation) -> Set

Gets the cocircular circumcenters of the Voronoi tessellation vorn as a Set of circumcenter indices. These are circumcenters that come from triangles that are cocircular with another adjoining triangle.

source
DelaunayTriangulation.get_generatorMethod
get_generator(vor::VoronoiTessellation, i) -> NTuple{2, Number}
get_generator(vor::VoronoiTessellation, i...) -> NTuple{length(i), NTuple{2, Number}}

Gets the coordinates for the generators i..., returned as Tuples of the form (x, y) for each generator.

source
DelaunayTriangulation.get_generatorsMethod
get_generators(vorn::VoronoiTessellation) -> Dict{Vertex,Point}

Gets the generators of the Voronoi tessellation vorn as a Dict, mapping vertices to their coordinates. These coordinates are given as Tuples of the form (x, y).

source
DelaunayTriangulation.get_polygonMethod
get_polygon(vor::VoronoiTessellation, i) -> Vector{Vertex}

Gets the vector of vertices corresponding to the ith polygon, given in counter-clockwise order and with the first and last vertices equal. To obtain the coordinates, see get_polygon_point.

source
DelaunayTriangulation.get_polygon_pointMethod
get_polygon_point(vor::VoronoiTessellation, i) -> NTuple{2, Number}
get_polygon_point(vor::VoronoiTessellation, i...) -> NTuple{length(i), NTuple{2, Number}}

Gets the coordinates corresponding to the vertices i... of the polygons, returned as Tuples of the form (x, y) for each vertex.

source
DelaunayTriangulation.get_polygon_pointsMethod
get_polygon_points(vorn::VoronoiTessellation) -> Vector{Point}

Gets the polygon points of the Voronoi tessellation vorn. These are the vertices of the Voronoi polygons, and are given as Tuples of the form (x, y).

source
DelaunayTriangulation.get_polygonsMethod
get_polygons(vorn::VoronoiTessellation) -> Dict{Index,Vector{Vertex}}

Gets the polygons of the Voronoi tessellation vorn as a Dict, mapping polygon indices to their vertices, where the vertices refer to points in get_polygon_points(vorn). The vertices are given in counter-clockwise order and the first and last vertices are equal.

source
DelaunayTriangulation.get_triangle_to_circumcenterMethod
get_triangle_to_circumcenter(vorn::VoronoiTessellation) -> Dict{Triangle,Index}

Gets the triangles of the Voronoi tessellation vorn as a Dict, mapping triangle indices to their corresponding circumcenters. The circumcenters are given as their vertices, referring to points in get_polygon_points(vorn).

source
DelaunayTriangulation.num_polygon_verticesMethod
num_polygon_vertices(vor::VoronoiTessellation) -> Integer

Returns the number of polygon vertices in the Voronoi tessellation vor. This might include duplicate vertices if get_polygon_points(vor) has duplicates.

source
DelaunayTriangulation.polygon_boundsFunction
polygon_bounds(vorn::VoronoiTessellation, unbounded_extension_factor=0.0; include_polygon_vertices=true) -> (Number, Number, Number, Number)

Gets the bounding box for the Voronoi tessellation vorn.

Arguments

  • vorn::VoronoiTessellation: The Voronoi tessellation.
  • unbounded_extension_factor=0.0: The factor by which to extend the bounding box for unbounded polygons.

Keyword Arguments

  • include_polygon_vertices=true: If true, then the bounding box will also include the polygon vertices. Otherwise, only the generators are included.

Output

  • xmin: Given by xmin′ - unbounded_extension_factor * (xmin′ - xmin′), where xmin′ is the original minimum x-coordinate of the computed bounding box and similarly for xmax′.
  • xmax: Given by xmax′ + unbounded_extension_factor * (xmax′ - xmax′), where xmax′ is the original maximum x-coordinate of the computed bounding box and similarly for xmin′.
  • ymin: Given by ymin′ - unbounded_extension_factor * (ymin′ - ymin′), where ymin′ is the original minimum y-coordinate of the computed bounding box and similarly for ymax′.
  • ymax: Given by ymax′ + unbounded_extension_factor * (ymax′ - ymax′), where ymax′ is the original maximum y-coordinate of the computed bounding box and similarly for ymin′.
source

Polygon

We sometimes find it useful to use a specific Polygon struct for representing polygons.

DelaunayTriangulation.PolygonType
Polygon{T,V,P} <: AbstractVector{T}

A struct for representing a polygon. The vertices are to be a counter-clockwise list of integers, where the integers themselves refer to points in points.

Fields

  • vertices::V: A list of integers that refer to points in points. The last vertex shokuld not be the same as the first.
  • points::P: A list of points.
Aliasing

In the case where vertices[begin] ≠ vertices[end], the vertices field is exactly the same as the input vertices. Where vertices[begin] = vertices[end], the vertices field is a view of vertices that excludes the last element.

source

ShuffledPolygonLinkedList

For computing constrained Delaunay triangulations, we use a ShuffledPolygonLinkedList struct to store the polygons involved.

DelaunayTriangulation.ShuffledPolygonLinkedListType
ShuffledPolygonLinkedList{I,T}

Data structure for representing a polygon as a doubly-linked list. In the descriptions below, π is used to denote the shuffled_indices vector.

Fields

  • next::Vector{I}: The next vertices, so that next[π[i]] is the vertex after S[π[i]].
  • prev::Vector{I}: The prev vertices, so that prev[π[i]] is the vertex before S[π[i]].
  • shuffled_indices::Vector{I}: The shuffled indices of the vertices, so that S[π[i]] is the ith vertex.
  • k::I: The number of vertices in the polygon.
  • S::T: The vertices of the polygon. This should not be a circular vector, i.e. S[begin] ≠ S[end], and must use one-based indexing. Additionally, the vertices must be provided in counter-clockwise order - this is NOT checked.

Constructor

To construct this, use

ShuffledPolygonLinkedList(S::Vector; rng::AbstractRNG=Random.default_rng())

The argument rng is used for shuffling the shuffled_indices vector.

source
DelaunayTriangulation.delete_vertex!Method
delete_vertex!(list::ShuffledPolygonLinkedList, i)

Deletes the vertex S[πᵢ] from the linked list, where πᵢ = list.shuffled_indices[i]. This deletion is done symbolically rather than by mutating the vectors in list. In particular, we perform

  • next[prev[πᵢ]] = next[πᵢ]
  • prev[next[πᵢ]] = prev[πᵢ]

which is the same as removing S[πᵢ] from the linked list.

source
DelaunayTriangulation.get_tripletMethod
get_triplet(list::ShuffledPolygonLinkedList, i) -> (Vertex, Vertex, Vertex)

Returns (u, v, w) = (S[πᵢ], S[next[πᵢ]], S[prev[πᵢ]]), where πᵢ = list.shuffled_indices[i] and S = list.S.

source
DelaunayTriangulation.reset!Method
reset!(list::ShuffledPolygonLinkedList; rng::AbstractRNG=Random.default_rng())

Resets the linked list, so that list.next[i] = mod1(i+1, list.k) and list.prev[i] = mod1(i-1, list.k), and also reshuffles the list.shuffled_indices vector.

source
DelaunayTriangulation.swap_permutation!Method
swap_permutation!(list::ShuffledPolygonLinkedList, i, j)

Reorders the permutation list.shuffled_indices of the linked list, swapping πᵢ and πⱼ where πₖ = list.shuffled_indices[k].

source

Points (Primitive Interface)

Here are functions that are used for defining and working with points in the package.

DelaunayTriangulation._getxMethod
_getx(p) -> Float64

Get the x-coordinate of p as a Float64.

Examples

julia> using DelaunayTriangulation

julia> p = (0.37, 0.7);

julia> DelaunayTriangulation._getx(p)
0.37

julia> p = (0.37f0, 0.7f0);

julia> DelaunayTriangulation._getx(p)
0.3700000047683716
source
DelaunayTriangulation._getxyMethod
_getxy(p) -> NTuple{2, Float64}

Get the coordinates of p as a Tuple.

Examples

julia> using DelaunayTriangulation

julia> p = [0.3, 0.5];

julia> DelaunayTriangulation._getxy(p)
(0.3, 0.5)

julia> p = [0.3f0, 0.5f0];

julia> DelaunayTriangulation._getxy(p)
(0.30000001192092896, 0.5)
source
DelaunayTriangulation._getyMethod
_gety(p) -> Float64

Get the y-coordinate of p as a Float64.

Examples

julia> using DelaunayTriangulation

julia> p = (0.5, 0.5);

julia> DelaunayTriangulation._gety(p)
0.5

julia> p = (0.5f0, 0.5f0);

julia> DelaunayTriangulation._gety(p)
0.5
source
DelaunayTriangulation.find_point_indexMethod
find_point_index(points, x, y) -> Integer
find_point_index(points, p) -> Integer

Returns an index of a point in points that is equal to p = (x, y). If no such point exists, then 0 is returned.

source
DelaunayTriangulation.getxMethod
getx(p) -> Number

Get the x-coordinate of p.

Examples

julia> using DelaunayTriangulation

julia> p = (0.3, 0.7);

julia> getx(p)
0.3
source
DelaunayTriangulation.getxyMethod
getxy(p) -> NTuple{2, Number}

Get the coordinates of p as a Tuple.

Examples

julia> using DelaunayTriangulation

julia> p = [0.9, 23.8];

julia> getxy(p)
(0.9, 23.8)
source
DelaunayTriangulation.getyMethod
gety(p) -> Number

Get the y-coordinate of p.

Examples

julia> using DelaunayTriangulation

julia> p = (0.9, 1.3);

julia> gety(p)
1.3
source
DelaunayTriangulation.lexicographic_orderMethod
lexicographic_order(points) -> Vector{Int}

Returns a set of indices that give the lexicographic ordering of points, meaning the indices so that the points are sorted first by x and then by y.

Examples

julia> using DelaunayTriangulation

julia> points = [(1.0, 5.0), (0.0, 17.0), (0.0, 13.0), (5.0, 17.3), (3.0, 1.0), (5.0, -2.0)]
6-element Vector{Tuple{Float64, Float64}}:
 (1.0, 5.0)
 (0.0, 17.0)
 (0.0, 13.0)
 (5.0, 17.3)
 (3.0, 1.0)
 (5.0, -2.0)

julia> DelaunayTriangulation.lexicographic_order(points)
6-element Vector{Int64}:
 3
 2
 1
 5
 6
 4

julia> hcat(points, points[ans])
6×2 Matrix{Tuple{Float64, Float64}}:
 (1.0, 5.0)   (0.0, 13.0)
 (0.0, 17.0)  (0.0, 17.0)
 (0.0, 13.0)  (1.0, 5.0)
 (5.0, 17.3)  (3.0, 1.0)
 (3.0, 1.0)   (5.0, -2.0)
 (5.0, -2.0)  (5.0, 17.3)
source
DelaunayTriangulation.mean_pointsFunction
mean_points(points[, vertices = each_point_index(points)]) -> NTuple{2, Number}

Returns the mean of the points in points indexed by vertices, given as a Tuple of the form (mean_x, mean_y).

Examples

julia> using DelaunayTriangulation

julia> points = [(1.0, 2.0), (2.3, 5.0), (17.3, 5.3)]
3-element Vector{Tuple{Float64, Float64}}:
 (1.0, 2.0)
 (2.3, 5.0)
 (17.3, 5.3)

julia> DelaunayTriangulation.mean_points(points)
(6.866666666666667, 4.1000000000000005)

julia> (1.0 + 2.3 + 17.3)/3, (2.0 + 5.0 + 5.3)/3
(6.866666666666667, 4.1000000000000005)

julia> points = [1.0 2.3 17.3; 2.0 5.0 5.3]
2×3 Matrix{Float64}:
 1.0  2.3  17.3
 2.0  5.0   5.3

julia> DelaunayTriangulation.mean_points(points)
(6.866666666666667, 4.1000000000000005)

julia> DelaunayTriangulation.mean_points(points, (1, 3))
(9.15, 3.65)

julia> (1.0 + 17.3)/2, (2.0 + 5.3)/2
(9.15, 3.65)
source
DelaunayTriangulation.points_are_uniqueMethod
points_are_unique(points) -> Bool

Returns true if all points in points are unique.

Examples

julia> using DelaunayTriangulation

julia> points = [1.0 2.0 3.0 4.0 5.0; 0.0 5.5 2.0 1.3 17.0]
2×5 Matrix{Float64}:
 1.0  2.0  3.0  4.0   5.0
 0.0  5.5  2.0  1.3  17.0

julia> DelaunayTriangulation.points_are_unique(points)
true

julia> points[:, 4] .= points[:, 1];

julia> DelaunayTriangulation.points_are_unique(points)
false
source

Edges (Primitive Interface)

Here are functions that are used for defining and working with edges in the package.

DelaunayTriangulation.add_edge!Method
add_edge!(E, e...)

Add the edges e... to E.

Examples

julia> using DelaunayTriangulation

julia> E = Set(((1,5),(17,10),(5,3)))
Set{Tuple{Int64, Int64}} with 3 elements:
  (5, 3)
  (17, 10)
  (1, 5)

julia> DelaunayTriangulation.add_edge!(E, (3, 2))

julia> E
Set{Tuple{Int64, Int64}} with 4 elements:
  (3, 2)
  (5, 3)
  (17, 10)
  (1, 5)

julia> DelaunayTriangulation.add_edge!(E, (1, -3), (5, 10), (1, -1))

julia> E
Set{Tuple{Int64, Int64}} with 7 elements:
  (3, 2)
  (5, 10)
  (1, -3)
  (1, -1)
  (5, 3)
  (17, 10)
  (1, 5)
source
DelaunayTriangulation.add_to_edges!Method
add_to_edges!(E, e)

Add the edge e to E.

Examples

julia> using DelaunayTriangulation

julia> E = Set(((1, 2),(3,5)))
Set{Tuple{Int64, Int64}} with 2 elements:
  (1, 2)
  (3, 5)

julia> DelaunayTriangulation.add_to_edges!(E, (1, 5))
Set{Tuple{Int64, Int64}} with 3 elements:
  (1, 2)
  (3, 5)
  (1, 5)
source
DelaunayTriangulation.compare_unoriented_edgesMethod
compare_unoriented_edges(u, v) -> Bool

Compare the unoriented edges u and v, i.e. compare the vertices of u and v in any order.

Examples

julia> using DelaunayTriangulation

julia> u = (1, 3);

julia> v = (5, 3);

julia> DelaunayTriangulation.compare_unoriented_edges(u, v)
false

julia> v = (1, 3);

julia> DelaunayTriangulation.compare_unoriented_edges(u, v)
true

julia> v = (3, 1);

julia> DelaunayTriangulation.compare_unoriented_edges(u, v)
true
source
DelaunayTriangulation.delete_edge!Method
delete_edge!(E, e...)

Delete the edges e... from E.

Examples

julia> using DelaunayTriangulation

julia> E = Set(([1,2],[10,15],[1,-1],[13,23],[1,5]))
Set{Vector{Int64}} with 5 elements:
  [10, 15]
  [1, 5]
  [1, 2]
  [1, -1]
  [13, 23]

julia> DelaunayTriangulation.delete_edge!(E, [10,15])

julia> E
Set{Vector{Int64}} with 4 elements:
  [1, 5]
  [1, 2]
  [1, -1]
  [13, 23]

julia> DelaunayTriangulation.delete_edge!(E, [1,5], [1, -1])

julia> E
Set{Vector{Int64}} with 2 elements:
  [1, 2]
  [13, 23]
source
DelaunayTriangulation.delete_from_edges!Method
delete_from_edges!(E, e)

Delete the edge e from E.

Examples

julia> using DelaunayTriangulation

julia> E = Set(([1,2],[5,15],[17,10],[5,-1]))
Set{Vector{Int64}} with 4 elements:
  [17, 10]
  [5, 15]
  [5, -1]
  [1, 2]

julia> DelaunayTriangulation.delete_from_edges!(E, [5, 15])
Set{Vector{Int64}} with 3 elements:
  [17, 10]
  [5, -1]
  [1, 2]
source
DelaunayTriangulation.edge_typeMethod
edge_type(E) -> DataType

Get the type of edges in E.

Examples

julia> using DelaunayTriangulation

julia> e = Set(((1,2),(2,3),(17,5)))
Set{Tuple{Int64, Int64}} with 3 elements:
  (1, 2)
  (17, 5)
  (2, 3)

julia> DelaunayTriangulation.edge_type(e)
Tuple{Int64, Int64}

julia> e = [[1,2],[3,4],[17,3]]
3-element Vector{Vector{Int64}}:
 [1, 2]
 [3, 4]
 [17, 3]

julia> DelaunayTriangulation.edge_type(e)
Vector{Int64} (alias for Array{Int64, 1})
source
DelaunayTriangulation.edge_verticesMethod
edge_vertices(e) -> NTuple{2, Vertex}

Get the vertices of e

Examples

julia> using DelaunayTriangulation

julia> e = (1, 5);

julia> edge_vertices(e)
(1, 5)

julia> e = [23, 50];

julia> edge_vertices(e)
(23, 50)
source
DelaunayTriangulation.initialMethod
initial(e) -> Vertex

Get the initial vertex of e.

Examples

julia> using DelaunayTriangulation

julia> e = (1, 3);

julia> DelaunayTriangulation.initial(e)
1

julia> e = [2, 5];

julia> DelaunayTriangulation.initial(e)
2
source
DelaunayTriangulation.num_edgesMethod
num_edges(E) -> Integer

Get the number of edges in E.

Examples

julia> using DelaunayTriangulation

julia> e = [(1, 2), (3, 4), (1, 5)];

julia> num_edges(e)
3
source
DelaunayTriangulation.reverse_edgeMethod
reverse_edge(e) -> Edge

Get the edge with the vertices of e in reverse order.

Examples

julia> using DelaunayTriangulation

julia> e = (17, 3);

julia> DelaunayTriangulation.reverse_edge(e)
(3, 17)

julia> e = [1, 2];

julia> DelaunayTriangulation.reverse_edge(e)
2-element Vector{Int64}:
 2
 1
source
DelaunayTriangulation.terminalMethod
terminal(e) -> Vertex

Get the terminal vertex of e.

Examples

julia> using DelaunayTriangulation

julia> e = (1, 7);

julia> DelaunayTriangulation.terminal(e)
7

julia> e = [2, 13];

julia> DelaunayTriangulation.terminal(e)
13
source

Triangles (Primitive Interface)

Here are functions that are used for defining and working with triangles in the package.

DelaunayTriangulation.add_to_triangles!Method
add_to_triangles!(T, V)

Add the triangle V to the collection of triangles T.

Examples

julia> using DelaunayTriangulation

julia> T = Set(((1, 2, 3), (17, 8, 9)));

julia> DelaunayTriangulation.add_to_triangles!(T, (1, 5, 12))
Set{Tuple{Int64, Int64, Int64}} with 3 elements:
  (1, 5, 12)
  (1, 2, 3)
  (17, 8, 9)

julia> DelaunayTriangulation.add_to_triangles!(T, (-1, 3, 6))
Set{Tuple{Int64, Int64, Int64}} with 4 elements:
  (1, 5, 12)
  (1, 2, 3)
  (17, 8, 9)
  (-1, 3, 6)
source
DelaunayTriangulation.compare_triangle_collectionsMethod
compare_triangle_collections(T, V) -> Bool

Compare the collections of triangles T and V by comparing their triangles according to compare_triangles.

Examples

julia> using DelaunayTriangulation

julia> T = Set(((1, 2, 3), (4, 5, 6), (7, 8, 9)));

julia> V = [[2, 3, 1], [4, 5, 6], [9, 7, 8]];

julia> DelaunayTriangulation.compare_triangle_collections(T, V)
true

julia> V[1] = [17, 19, 20];

julia> DelaunayTriangulation.compare_triangle_collections(T, V)
false

julia> V = [[1, 2, 3], [8, 9, 7]];

julia> DelaunayTriangulation.compare_triangle_collections(T, V)
false
source
DelaunayTriangulation.compare_trianglesMethod
compare_triangles(T, V) -> Bool

Compare the triangles T and V by comparing their vertices up to rotation.

Examples

julia> using DelaunayTriangulation

julia> T1 = (1, 5, 10);

julia> T2 = (17, 23, 20);

julia> DelaunayTriangulation.compare_triangles(T1, T2)
false

julia> T2 = (5, 10, 1);

julia> DelaunayTriangulation.compare_triangles(T1, T2)
true

julia> T2 = (10, 1, 5);

julia> DelaunayTriangulation.compare_triangles(T1, T2)
true

julia> T2 = (10, 5, 1);

julia> DelaunayTriangulation.compare_triangles(T1, T2)
false
source
DelaunayTriangulation.construct_positively_oriented_triangleMethod
construct_positively_oriented_triangle(::Type{V}, i, j, k, points) where {V} -> V

Construct a triangle of type V from vertices i, j, and k such that the triangle is positively oriented, using points for the coordinates.

Examples

julia> using DelaunayTriangulation

julia> points = [(0.0, 0.0), (0.0, 1.0), (1.0, 0.0)];

julia> DelaunayTriangulation.construct_positively_oriented_triangle(NTuple{3, Int}, 1, 2, 3, points)
(2, 1, 3)

julia> DelaunayTriangulation.construct_positively_oriented_triangle(NTuple{3, Int}, 2, 3, 1, points)
(3, 2, 1)

julia> DelaunayTriangulation.construct_positively_oriented_triangle(NTuple{3, Int}, 2, 1, 3, points)
(2, 1, 3)

julia> DelaunayTriangulation.construct_positively_oriented_triangle(NTuple{3, Int}, 3, 2, 1, points)
(3, 2, 1)

julia> points = [(1.0, 1.0), (2.5, 2.3), (17.5, 23.0), (50.3, 0.0), (-1.0, 2.0), (0.0, 0.0), (5.0, 13.33)];

julia> DelaunayTriangulation.construct_positively_oriented_triangle(Vector{Int}, 5, 3, 2, points)
3-element Vector{Int64}:
 3
 5
 2

julia> DelaunayTriangulation.construct_positively_oriented_triangle(Vector{Int}, 7, 1, 2, points)
3-element Vector{Int64}:
 7
 1
 2

julia> DelaunayTriangulation.construct_positively_oriented_triangle(Vector{Int}, 7, 2, 1, points)
3-element Vector{Int64}:
 2
 7
 1

julia> DelaunayTriangulation.construct_positively_oriented_triangle(Vector{Int}, 5, 4, 3, points)
3-element Vector{Int64}:
 5
 4
 3
source
DelaunayTriangulation.delete_from_triangles!Method
delete_from_triangles!(T, V)

Delete the triangle V from the collection of triangles T. Only deletes V if V is in T up to rotation.

Examples

julia> using DelaunayTriangulation

julia> V = Set(((1, 2, 3), (4, 5, 6), (7, 8, 9)))
Set{Tuple{Int64, Int64, Int64}} with 3 elements:
  (7, 8, 9)
  (4, 5, 6)
  (1, 2, 3)

julia> DelaunayTriangulation.delete_from_triangles!(V, (4, 5, 6))
Set{Tuple{Int64, Int64, Int64}} with 2 elements:
  (7, 8, 9)
  (1, 2, 3)

julia> DelaunayTriangulation.delete_from_triangles!(V, (9, 7, 8))
Set{Tuple{Int64, Int64, Int64}} with 1 element:
  (1, 2, 3)
source
DelaunayTriangulation.getiMethod
geti(T) -> Vertex

Get the first vertex of T.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.geti((1, 2, 3))
1

julia> DelaunayTriangulation.geti([2, 5, 1])
2
source
DelaunayTriangulation.getjMethod
getj(T) -> Vertex

Get the second vertex of T.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.getj((5, 6, 13))
6

julia> DelaunayTriangulation.getj([10, 19, 21])
19
source
DelaunayTriangulation.getkMethod
getk(T) -> Vertex

Get the third vertex of T.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.getk((1,2,3))
3

julia> DelaunayTriangulation.getk([1,2,3])
3
source
DelaunayTriangulation.num_trianglesMethod
num_triangles(T) -> Integer

Get the number of triangles in T.

Examples

julia> using DelaunayTriangulation

julia> T1, T2, T3 = (1, 5, 10), (17, 23, 10), (-1, 10, 5);

julia> T = Set((T1, T2, T3));

julia> num_triangles(T)
3
source
DelaunayTriangulation.rotate_triangleMethod
rotate_triangle(T, rotation) -> Triangle

Rotate the vertices of T by rotation. In particular, if T = (i, j, k):

  • rotation = 0: (i, j, k)
  • rotation = 1: (j, k, i)
  • rotation = 2: (k, i, j)
  • Otherwise, return rotate_triangle(T, rotation % 3).

Examples

julia> using DelaunayTriangulation

julia> T = (1, 2, 3)
(1, 2, 3)

julia> DelaunayTriangulation.rotate_triangle(T, 0)
(1, 2, 3)

julia> DelaunayTriangulation.rotate_triangle(T, 1)
(2, 3, 1)

julia> DelaunayTriangulation.rotate_triangle(T, 2)
(3, 1, 2)

julia> DelaunayTriangulation.rotate_triangle(T, 3)
(1, 2, 3)
source
DelaunayTriangulation.sort_triangleFunction
sort_triangle(T) -> Triangle
sort_triangle(i, j, k) -> Triangle

Sort the triangle T = (i, j, k) so that its last vertex is the smallest, respecting the orientation of T.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.sort_triangle((1, 5, 3))
(5, 3, 1)

julia> DelaunayTriangulation.sort_triangle((1, -1, 2))
(2, 1, -1)

julia> DelaunayTriangulation.sort_triangle((3, 2, 1))
(3, 2, 1)
source
DelaunayTriangulation.sort_trianglesMethod
sort_triangles(T) -> Triangle

Sort the triangles in T so that the first vertex of each triangle is the largest, respecting the orientation of the triangles. See sort_triangle.

Examples

julia> using DelaunayTriangulation

julia> T = Set(((1, 3, 2), (5, 2, 3), (10, 1, 13), (-1, 10, 12), (10, 1, 17), (5, 8, 2)))
Set{Tuple{Int64, Int64, Int64}} with 6 elements:
  (5, 8, 2)
  (10, 1, 13)
  (10, 1, 17)
  (5, 2, 3)
  (1, 3, 2)
  (-1, 10, 12)

julia> DelaunayTriangulation.sort_triangles(T)
Set{Tuple{Int64, Int64, Int64}} with 6 elements:
  (13, 10, 1)
  (3, 5, 2)
  (10, 12, -1)
  (5, 8, 2)
  (17, 10, 1)
  (3, 2, 1)
source
DelaunayTriangulation.triangle_typeMethod
triangle_type(::Type{T}) -> DataType

Get the triangle type of T.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.triangle_type(Set{NTuple{3,Int64}})
Tuple{Int64, Int64, Int64}

julia> DelaunayTriangulation.triangle_type(Vector{NTuple{3,Int32}})
Tuple{Int32, Int32, Int32}

julia> DelaunayTriangulation.triangle_type(Vector{Vector{Int64}})
Vector{Int64} (alias for Array{Int64, 1})
source
DelaunayTriangulation.triangle_verticesMethod
triangle_vertices(T) -> NTuple{3, Vertex}

Returns the vertices of T as a Tuple.

Examples

julia> using DelaunayTriangulation

julia> triangle_vertices((1, 5, 17))
(1, 5, 17)

julia> triangle_vertices([5, 18, 23]) # -> tuple
(5, 18, 23)
source

Boundary Nodes (Primitive Interface)

Here are functions that are used for defining and working with boundary nodes in the package.

DelaunayTriangulation.construct_boundary_edge_mapMethod
construct_boundary_edge_map(boundary_nodes::A, IntegerType::Type{I}=number_type(boundary_nodes), EdgeType::Type{E}=NTuple{2,IntegerType}) where {A,I,E} -> Dict

Constructs a map that takes boundary edges in boundary_nodes to a Tuple giving the edge's position in boundary_nodes. In particular, if dict = construct_boundary_edge_map(boundary_nodes), then dict[e] = (pos, ℓ) so that bn = get_boundary_nodes(boundary_nodes, pos) gives the boundary nodes associated with the section that e lives on, and get_boundary_nodes(bn, ℓ) is the first vertex of e.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.construct_boundary_edge_map([17, 18, 15, 4, 3, 17])
Dict{Tuple{Int64, Int64}, Tuple{Vector{Int64}, Int64}} with 5 entries:
  (18, 15) => ([17, 18, 15, 4, 3, 17], 2)
  (3, 17)  => ([17, 18, 15, 4, 3, 17], 5)
  (17, 18) => ([17, 18, 15, 4, 3, 17], 1)
  (4, 3)   => ([17, 18, 15, 4, 3, 17], 4)
  (15, 4)  => ([17, 18, 15, 4, 3, 17], 3)

julia> DelaunayTriangulation.construct_boundary_edge_map([[5, 17, 3, 9], [9, 18, 13, 1], [1, 93, 57, 5]])
Dict{Tuple{Int64, Int64}, Tuple{Int64, Int64}} with 9 entries:
  (18, 13) => (2, 2)
  (17, 3)  => (1, 2)
  (9, 18)  => (2, 1)
  (13, 1)  => (2, 3)
  (3, 9)   => (1, 3)
  (93, 57) => (3, 2)
  (5, 17)  => (1, 1)
  (57, 5)  => (3, 3)
  (1, 93)  => (3, 1)

julia> DelaunayTriangulation.construct_boundary_edge_map([[[2, 5, 10], [10, 11, 2]], [[27, 28, 29, 30], [30, 31, 85, 91], [91, 92, 27]]])
Dict{Tuple{Int64, Int64}, Tuple{Tuple{Int64, Int64}, Int64}} with 12 entries:
  (92, 27) => ((2, 3), 2)
  (2, 5)   => ((1, 1), 1)
  (11, 2)  => ((1, 2), 2)
  (10, 11) => ((1, 2), 1)
  (30, 31) => ((2, 2), 1)
  (91, 92) => ((2, 3), 1)
  (29, 30) => ((2, 1), 3)
  (31, 85) => ((2, 2), 2)
  (27, 28) => ((2, 1), 1)
  (5, 10)  => ((1, 1), 2)
  (28, 29) => ((2, 1), 2)
  (85, 91) => ((2, 2), 3)
source
DelaunayTriangulation.construct_ghost_vertex_mapMethod
construct_ghost_vertex_map(boundary_nodes::A, IntegerType::Type{I}=number_type(boundary_nodes)) where {A,I} -> Dict

Given a set of boundary_nodes, returns a Dict that maps ghost vertices to their associated section in boundary_nodes. There are three cases:

  • has_multiple_curves(boundary_nodes)

Returns dict::Dict{I, NTuple{2, I}}, mapping ghost vertices i to Tuples (m, n) so that get_boundary_nodes(boundary_nodes, m, n) are the boundary nodes associated with i, i.e. the nth section of the mth curve is associated with the ghost vertex i.

  • has_multiple_sections(boundary_nodes)

Returns dict::Dict{I, I}, mapping ghost vertices i to n so that get_boundary_nodes(boundary_nodes, n) are the boundary nodes associated with i, i.e. the nth section of the boundary is associated with the ghost vertex i.

  • otherwise

Returns dict::Dict{I, A}, mapping the ghost vertex i to boundary_nodes.

Examples

julia> using DelaunayTriangulation

julia> gv_map = DelaunayTriangulation.construct_ghost_vertex_map([1, 2, 3, 4, 5, 1])
Dict{Int64, Vector{Int64}} with 1 entry:
  -1 => [1, 2, 3, 4, 5, 1]

julia> gv_map = DelaunayTriangulation.construct_ghost_vertex_map([[17, 29, 23, 5, 2, 1], [1, 50, 51, 52], [52, 1]])
Dict{Int64, Int64} with 3 entries:
  -1 => 1
  -3 => 3
  -2 => 2

julia> gv_map = DelaunayTriangulation.construct_ghost_vertex_map([[[1, 5, 17, 18, 1]], [[23, 29, 31, 33], [33, 107, 101], [101, 99, 85, 23]]])
Dict{Int64, Tuple{Int64, Int64}} with 4 entries:
  -1 => (1, 1)
  -3 => (2, 2)
  -2 => (2, 1)
  -4 => (2, 3)

Extended help

This map can be useful for iterating over all boundary nodes. For example, you can iterate over all sections of a boundary using:

gv_map = construct_ghost_vertex_map(boundary_nodes)
for (ghost_vertex, section) in gv_map 
    nodes = get_boundary_nodes(boundary_nodes, section)
    # do something with nodes
end

This works for any form of boundary_nodes.

source
DelaunayTriangulation.construct_ghost_vertex_rangesMethod
construct_ghost_vertex_ranges(boundary_nodes::A, IntegerType::Type{I}=number_type(boundary_nodes)) where {A,I} -> Dict

Given a set of boundary_nodes, returns a Dict that maps ghost vertices to the range of all ghost vertices that the corresponding boundary curve could correspond to.

Examples

julia> using DelaunayTriangulation

julia> boundary_nodes = [
                  [
                      [1, 2, 3, 4], [4, 5, 6, 1]
                  ],
                  [
                      [18, 19, 20, 25, 26, 30]
                  ],
                  [
                      [50, 51, 52, 53, 54, 55], [55, 56, 57, 58], [58, 101, 103, 105, 107, 120], [120, 121, 122, 50]
                  ]
              ]
3-element Vector{Vector{Vector{Int64}}}:
 [[1, 2, 3, 4], [4, 5, 6, 1]]
 [[18, 19, 20, 25, 26, 30]]
 [[50, 51, 52, 53, 54, 55], [55, 56, 57, 58], [58, 101, 103, 105, 107, 120], [120, 121, 122, 50]]

julia> DelaunayTriangulation.construct_ghost_vertex_ranges(boundary_nodes)
Dict{Int64, UnitRange{Int64}} with 7 entries:
  -5 => -7:-4
  -1 => -2:-1
  -7 => -7:-4
  -3 => -3:-3
  -2 => -2:-1
  -4 => -7:-4
  -6 => -7:-4
source
DelaunayTriangulation.delete_boundary_node!Method
delete_boundary_node!(boundary_nodes, pos)

Deletes the boundary node at the position pos from boundary_nodes. Here, pos[1] is such that get_boundary_nodes(boundary_nodes, pos[1]) is the section that the node will be deleted from, and pos[2] gives the position of the array to delete from. In particular,

delete_boundary_node!(boundary_nodes, pos)

is the same as

deleteat!(get_boundary_nodes(boundary_nodes, pos[1]), pos[2])

Examples

julia> using DelaunayTriangulation

julia> boundary_nodes = [71, 25, 33, 44, 55, 10];

julia> DelaunayTriangulation.delete_boundary_node!(boundary_nodes, (boundary_nodes, 4))
5-element Vector{Int64}:
 71
 25
 33
 55
 10

julia> boundary_nodes = [[7, 13, 9, 25], [25, 26, 29, 7]];

julia> DelaunayTriangulation.delete_boundary_node!(boundary_nodes, (2, 3))
2-element Vector{Vector{Int64}}:
 [7, 13, 9, 25]
 [25, 26, 7]

julia> boundary_nodes = [[[17, 23, 18, 25], [25, 26, 81, 91], [91, 101, 17]], [[1, 5, 9, 13], [13, 15, 1]]];

julia> DelaunayTriangulation.delete_boundary_node!(boundary_nodes, ((2, 2), 2))
2-element Vector{Vector{Vector{Int64}}}:
 [[17, 23, 18, 25], [25, 26, 81, 91], [91, 101, 17]]
 [[1, 5, 9, 13], [13, 1]]
source
DelaunayTriangulation.each_boundary_nodeMethod
each_boundary_node(boundary_nodes) -> Iterator

Returns an iterator that goves over each node in boundary_nodes. It is assumed that boundary_nodes represents a single section or a single contiguous boundary; if you do want to loop over every boundary nodes for a boundary with multiple sections, you should to see the result from construct_ghost_vertex_map.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.each_boundary_node([7, 8, 19, 2, 17])
5-element Vector{Int64}:
  7
  8
 19
  2
 17

julia> DelaunayTriangulation.each_boundary_node([7, 8, 19, 2, 17, 7])
6-element Vector{Int64}:
  7
  8
 19
  2
 17
  7
source
DelaunayTriangulation.get_skeletonMethod
get_skeleton(boundary_nodes, IntegerType) -> empty(boundary_nodes)

Given a set of boundary nodes, returns the empty skeleton of the boundary nodes. This is essentially empty applied to boundary_nodes, but with vertices of type IntegerType. This is mainly needed for convert_boundary_curves!. You will need to implement a new method for this if you want your custom boundary node interface to be supported for curve-bounded domains.

source
DelaunayTriangulation.insert_boundary_node!Method
insert_boundary_node!(boundary_nodes, pos, node)

Inserts a boundary node node into boundary_nodes at the position pos. Here, pos[1] is such that get_boundary_nodes(boundary_nodes, pos[1]) is the section that the node will be inserted onto, and pos[2] gives the position of the array to insert node into. In particular,

insert_boundary_node!(boundary_nodes, pos, node)

is the same as

insert!(get_boundary_nodes(boundary_nodes, pos[1]), pos[2], node)

Examples

julia> using DelaunayTriangulation

julia> boundary_nodes = [1, 2, 3, 4, 5, 1];

julia> DelaunayTriangulation.insert_boundary_node!(boundary_nodes, (boundary_nodes, 4), 23)
7-element Vector{Int64}:
  1
  2
  3
 23
  4
  5
  1

julia> boundary_nodes = [[7, 13, 9, 25], [25, 26, 29, 7]];

julia> DelaunayTriangulation.insert_boundary_node!(boundary_nodes, (2, 1), 57)
2-element Vector{Vector{Int64}}:
 [7, 13, 9, 25]
 [57, 25, 26, 29, 7]

julia> boundary_nodes = [[[17, 23, 18, 25], [25, 26, 81, 91], [91, 101, 17]], [[1, 5, 9, 13], [13, 15, 1]]];

julia> DelaunayTriangulation.insert_boundary_node!(boundary_nodes, ((1, 3), 3), 1001)
2-element Vector{Vector{Vector{Int64}}}:
 [[17, 23, 18, 25], [25, 26, 81, 91], [91, 101, 1001, 17]]
 [[1, 5, 9, 13], [13, 15, 1]]
source
DelaunayTriangulation.num_boundary_edgesMethod
num_boundary_edges(boundary_nodes) -> Integer

Get the number of boundary edges in boundary_nodes, assuming that boundary_nodes defines a boundary with only one curve and a single section.

source
DelaunayTriangulation.num_curvesMethod
num_curves(boundary_nodes) -> Integer

Get the number of curves in boundary_nodes.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.num_curves([1, 2, 3, 1])
1

julia> DelaunayTriangulation.num_curves([[1, 2, 3], [3, 4, 1]])
1

julia> DelaunayTriangulation.num_curves([[[1, 2, 3], [3, 4, 1]], [[5, 6, 7, 8, 5]]])
2
source
DelaunayTriangulation.num_sectionsMethod
num_sections(boundary_nodes) -> Integer

Assuming boundary_nodes has only one curve, get the number of sections in boundary_nodes.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.num_sections([1, 2, 3, 4, 5, 1])
1

julia> DelaunayTriangulation.num_sections([[1, 2, 3, 4], [4, 5, 1]])
2

julia> DelaunayTriangulation.num_sections([[1, 2, 3], [3, 4, 5, 6, 7, 8], [8, 9], [9, 1]])
4
source