Data Structures
DelaunayTriangulation.Adjacent
— TypeAdjacent{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})
DelaunayTriangulation.Adjacent2Vertex
— TypeAdjacent2Vertex{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})
DelaunayTriangulation.Graph
— TypeGraph{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}})
DelaunayTriangulation.Triangulation
— TypeTriangulation{P,T,BN,W,I,E,Es,BC,BCT,BEM,GVM,GVR,BPL,C,BE}
Struct representing a triangulation, as constructed by triangulate
.
Accessing the fields themselves using e.g. tri.field
is not recommended and is not intended to be in the public API. You should be using the accessor functions, e.g. instead of tri.points
do get_points(tri)
. Similarly, for the iterators, e.g. tri.triangles
, each_triangle(tri)
is recommended instead.
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 with 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 i
th vertex. The weights should have the same type as the coordinates in points
.
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 i
th element of the Tuple
is associated with the ghost vertex -i
, i.e. the i
th section as indicated by ghost_vertex_map
. If the i
th 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!
.
DelaunayTriangulation.VoronoiTessellation
— TypeVoronoiTessellation{Tr<:Triangulation,P,I,T,S,E}
Struct for representing a Voronoi tessellation.
See also voronoi
.
Accessing the fields themselves using e.g. vorn.field
is not recommended and is not intended to be in the public API. You should be using the accessor functions, e.g. instead of vorn.adjacent
do get_adjacent(vorn)
. Similarly, for the iterators, e.g. vorn.generators
, each_generators(vorn)
is recommended instead.
In the case that the underlying triangulation is weighted, then this struct represents the power diagram, and instead of circumcenters the points are orthocenters computed with triangle_orthocenter
.
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}
: ADict
that maps vertices of generators to coordinates. These are simply the points present in the triangulation. ADict
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 alsoget_polygon_coordinates
.)polygons::Dict{I,Vector{I}}
: ADict
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}
: ADict
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}
: ADict
mapping a triangle to its circumcenter index. The triangles are sorted such that the minimum vertex is last.unbounded_polygons::Set{I}
: ASet
of indices of the unbounded polygons.cocircular_circumcenters::S
: ASet
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}
: ASet
of indices of the polygons that are on the boundary of the tessellation. Only relevant for clipped tessellations, otherwise seeunbounded_polygons
.
DelaunayTriangulation.ConvexHull
— TypeConvexHull{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 thatvertices[begin] == vertices[end]
.
Constructors
ConvexHull(points, vertices)
convex_hull(points; predicates=AdaptiveKernel(), IntegerType=Int)
DelaunayTriangulation.InsertionEventHistory
— TypeInsertionEventHistory{T,E}
A data structure for storing the changes to the triangulation during the insertion of a point.
Fields
added_triangles::Set{T}
: The triangles that were added.deleted_triangles::Set{T}
: The triangles that were deleted.added_segments::Set{E}
: The interior segments that were added.deleted_segments::Set{E}
: The interior segments that were deleted.added_boundary_segments::Set{E}
: The boundary segments that were added.deleted_boundary_segments::Set{E}
: The boundary segments that were deleted.
Constructor
The default constructor is available, but we also provide
InsertionEventHistory(tri::Triangulation)
which will initialise this struct with empty, appropriately sizehint!
ed, sets.
DelaunayTriangulation.PointLocationHistory
— TypePointLocationHistory{T,E,I}
History from using find_triangle
.
Fields
triangles::Vector{T}
: The visited triangles.collinear_segments::Vector{E}
: Segments collinear with the original linepq
using to jump.collinear_point_indices::Vector{I}
: This field contains indices to segments incollinear_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 ofpq
.right_verices::Vector{I}
: Vertices from the visited triangles to the right ofpq
.