Triangulation

The most important structure of the package is the Triangulation data structure. Its complete definition is given below, and then afterwards we give all the docstrings for methods that are useful for working with triangulations.

DelaunayTriangulation.TriangulationType
Triangulation{P,Ts,I,E,Es,BN,BNM,B,BIR,BPL}

Struct representing a Delaunay triangulation, as constructed via e.g. triangulate or generate_mesh.

Fields

  • points::P

The points used to construct the triangulation. Do note that this is not necessarly the same as the set of points appearing in the triangulation, as could happen e.g. if some points were deleted from the triangulation. If that is what you need, then see each_vertex.

  • triangles::Ts

The triangles in the triangulation, all given with positive orientation. This field can include ghost triangles. If you don't want ghost triangles, see each_solid_triangle.

  • adjacent::Adjacent{I, E}

The Adjacent map.

  • adjacent2vertex::Adjacent2Vertex{I,Es,E}

The Adjacent2Vertex map.

  • graph::Graph{I}

The Graph.

  • boundary_nodes::BN

The boundary nodes in the triangulation that you provide, matching the specification given in the documentation.

  • boundary_edge_map::BNM

This is a Dict that maps all of the boundary edges to their position in boundary_nodes. See also construct_boundary_edge_map.

  • boundary_index_ranges::BIR

This is an OrderedDict that maps a boundary index to a range of all other boundary indices that the corresponding boundary curve could correspond to. For example, for a curve with four segments, there are four possible boundary indices that a segment could correspond to. See also construct_boundary_index_ranges.

  • constrained_edges::Es

This is the set of extra edges added into the triangulation that you have provided. This will not include any of the other constrained edges from boundary_nodes.

  • all_constrained_edges::Es

This is the set of all constrained edges appearing in the triangulation, essentially given as the union of constrained_edges and boundary_nodes. This is different from the constrained_edges field as it contains edges from both of these fields, e.g. there might be an edge in constrained_edges that is not yet in the triangulation, but it will definitely not appear in all_constrained_edges.

  • convex_hull::ConvexHull{P,Vector{I}}

The ConvexHull of points. See also convex_hull.

  • representative_point_list::RPL

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

Constructors

There are several ways to construct this struct directly, although in most cases you should be using triangulate or generate_mesh.

Default Constructor

The default constructor is available, i.e.

Triangulation(
    points,
    triangles,
    adjacent,
    adjacent2vertex,
    graph,
    boundary_nodes,
    boundary_edge_map,
    boundary_map,
    boundary_index_ranges,
    constrained_edges,
    all_constrained_edges,
    convex_hull,
    representative_point_list
)

Empty Triangulation

An empty triangulation can be initalised with the following method,

Triangulation(points;
    IntegerType=Int,
    EdgeType=NTuple{2,IntegerType},
    TriangleType=NTuple{3,IntegerType},
    EdgesType=Set{EdgeType},
    TrianglesType=Set{TriangleType},
    boundary_nodes=IntegerType[],
    constrained_edges=initialise_edges(EdgesType),
    representative_point_list=get_empty_representative_points(IntegerType, number_type(points))
)

Triangulation From an Existing Mesh

A method is available from constructing a mesh from an existing set of points, triangles, and boundary nodes, mainly existing for the purpose of generate_mesh:

Triangulation(points, triangles, boundary_nodes;
    IntegerType=Int,
    EdgeType=NTuple{2,IntegerType},
    TriangleType=NTuple{3,IntegerType},
    EdgesType=Set{EdgeType},
    TrianglesType=Set{TriangleType},
    add_ghost_triangles=false
)

with add_ghost_triangles calling add_ghost_triangles! at the end of the constructor.

Wrapping a Triangulation with constrainededgepoints

This method is used in triangulate for wrapping a triangulation with a set of constraints:

remake_triangulation_with_constraints(triangulation, edges, boundary_nodes)

which returns (bn_map, bn_range, tri), where bn_map is the boundary_map that isn't yet added, bn_range is the boundary_index_range that isn't yet added either, and tri is the wrapped triangulation that includes the constrained edges and the boundary_nodes. Note that either edges or boundary_nodes can be nothing.

This can be used together with the method

replace_boundary_dict_information(triangulation, bn_map, bn_range)

which returns a new triangulation with the boundary_map and boundary_index_range replaced, completing the wrapper.

Each field has its own accessor:

There are several useful methods available for working with triangulations. We split these into several sections.

Adjacent Methods

DelaunayTriangulation.get_adjacentMethod
get_adjacent(tri, uv; check_existence=Val(has_multiple_segments(tri)))
get_adjacent(tri, u, v; check_existence=Val(has_multiple_segments(tri)))

Returns get_adjacent(tri, uv; check_existence) or get_adjacent(tri, u, v; check_existence).

DelaunayTriangulation.add_adjacent!Method
add_adjacent!(tri::Triangulation, uv, w)
add_adjacent!(tri::Triangulation, u, v, w)

Calls add_adjacent!(get_adjacent(tri), uv, w) or add_adjacent!(get_adjacent(tri), u, v, w), adding the edge uv = (u, v) to the adjacent map of tri with corresponding vertex w.

DelaunayTriangulation.delete_adjacent!Method
delete_adjacent!(tri::Triangulation, uv)
delete_adjacent!(tri::Triangulation, u, v)

Calls delete_adjacent!(get_adjacent(tri), uv) or delete_adjacent!(get_adjacent(tri), u, v), removing the edge uv = (u, v) from the adjacent map.

Adjacent2Vertex Methods

DelaunayTriangulation.add_adjacent2vertex!Method
add_adjacent2vertex!(tri::Triangulation, w, uv)
add_adjacent2vertex!(tri::Triangulation, w, u, v)

Calls add_adjacent2vertex!(get_adjacent2vertex(tri), w, uv) or add_adjacent2vertex!(get_adjacent2vertex(tri), w, u, v), pushing the edge uv = (u, v) into get_adjacent2vertex(tri, w).

DelaunayTriangulation.delete_adjacent2vertex!Method
delete_adjacent2vertex!(tri::Triangulation, w, uv)
delete_adjacent2vertex!(tri::Triangulation, w, u, v)

Calls delete_adjacent2vertex!(get_adjacent2vertex(tri), w, uv) or delete_adjacent2vertex!(get_adjacent2vertex(tri), w, u, v), deleting the edge uv = (u, v) from get_adjacent2vertex(tri, w).

Boundary Nodes Methods

DelaunayTriangulation.get_boundary_nodesMethod
get_boundary_nodes(bn, mnℓ...)

Get the boundary nodes from bn corresponding to the specified indices. See getboundarynodes.

get_boundary_nodes(tri::Triangulation, mnℓ...)

Returns get_boundary_nodes(get_boundary_nodes(tri), mnℓ...).

DelaunayTriangulation.map_boundary_indexMethod
map_boundary_index(tri::Triangulation, i)

Returns map_boundary_index(get_boundary_map(tri), i), mapping the boundary index i to the corresponding index in the boundary nodes.

DelaunayTriangulation.get_curve_indexMethod
get_curve_index(tri::Triangulation, i)

Returns get_curve_index(get_boundary_map(tri), i), the curve corresponding to the boundary index i.

DelaunayTriangulation.get_right_boundary_nodeMethod
get_right_boundary_node(tri::Triangulation, k, boundary_index)

Returns get_right_boundary_node(get_adjacent(tri), k, boundary_index, get_boundary_index_ranges(tri), Val(has_multiple_segments(tri))), the node to the right of the boundary node k. It is assumed that k is on the part of the boundary of tri with index boundary_index.

DelaunayTriangulation.get_left_boundary_nodeMethod
get_left_boundary_node(tri::Triangulation, k, boundary_index)

Returns get_left_boundary_node(get_adjacent(tri), k, boundary_index, get_boundary_index_ranges(tri), Val(has_multiple_segments(tri))), the node to the left of the boundary node k. It is assumed that k is on the part of the boundary of tri with index boundary_index.

DelaunayTriangulation.get_boundary_index_rangeMethod
get_boundary_index_range(tri::Triangulation, i)

Returns map_boundary_index(get_boundary_index_ranges(tri), i), the list of boundary indices belonging to the curve corresponding to the boundary index i.

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

Returns get_boundary_edge_map(get_boundary_nodes(tri), ij) or get_boundary_edge_map(get_boundary_nodes(tri), construct_edge(edge_type(tri), i, j)), returning a Tuple pos such that get_boundary_nodes(get_boundary_nodes(tri, pos[1]), pos[2]) is the edge ij = (i, j).

DelaunayTriangulation.split_boundary_edge!Method
split_boundary_edge!(tri::Triangulation, edge, node)
split_boundary_edge!(tri::Triangulation, i, j, node)

Given a triangulation tri, an edge = (i, j), and a node on edge, this function splits the boundary edge into two edges edge = (i, node) and (node, j), updating tri.boundary_nodes and tri.boundary_edge_map accordingly.

See also split_boundary_edge_at_collinear_segments!.

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

Given a triangulation tri and a list of collinear segments collinear_segments, assumed to represent a boundary edge (initial(first(collinear_segments)), terminal(last(collinear_segments))), splits the edge at the collinear_segments and updates tri.boundary_nodes and tri.boundary_edge_map accordingly.

See also split_boundary_edge!.

DelaunayTriangulation.contains_boundary_edgeMethod
contains_boundary_edge(tri::Triangulation, e)
contains_boundary_edge(tri::Triangulation, i, j)

Tests if the triangulation tri has the constrained boundary edge e = (i, j), returning e ∈ keys(get_boundary_edge_map(tri)).

DelaunayTriangulation.merge_constrained_edgesMethod
merge_constrained_edges(bn_map, boundary_nodes, constrained_edges::Es)
merge_constrained_edges(tri::Triangulation, bn_map=get_boundary_map(tri))

Merges the boundary edges, defined via boundary_nodes or get_boundary_nodes(tri), and the constrained_edges into a single collection. bn_map is used to iterate over all the boundary edges.

Missing docstring.

Missing docstring for merge_boundary_node!(::Triangulation, ::Any, ::Any). Check Documenter's build log for details.

Convex Hull Methods

DelaunayTriangulation.convex_hull!Method
convex_hull!(tri::Triangulation; reconstruct=is_constrained(tri))

Computes the convex hull of the points included in the triangulation tri. If reconstruct, then the method in convex_hull will be used, whereas if !reconstruct then the triangulation's boundary is extracted to get the convex hull directly. Note that this latter method fails if there are any constrained edges on the boundary of tri.

Edges Methods

DelaunayTriangulation.edge_typeMethod
edge_type(tri::Triangulation)

Returns the type used for representing edges in tri.

edge_type(vor::VoronoiTessellation)

Returns the type of the edges in the Voronoi tessellation vor.

DelaunayTriangulation.each_constrained_edgeMethod
each_constrained_edge(tri::Triangulation)

Returns each_edge(get_all_constrained_edges(tri)), an iterator over all constrained_edges (including all constrained boundary edges)`.

DelaunayTriangulation.contains_constrained_edgeMethod
contains_constrained_edge(tri::Triangulation, e)
contains_constrained_edge(tri::Triangulation, i, j)

Tests if tri contains the constrained edge e = (i, j) (or reverse_edge(e) = (j, i)).

DelaunayTriangulation.split_constrained_edge!Method
split_constrained_edge!(tri::Triangulation, constrained_edge, collinear_segments)

Calls split_constrained_edge!(get_constrained_edges(tri), constrained_edge, collinear_segments), splitting the constrained_edge which is assumed to represent the union of the collinear_segments so that we instead store those collinear_segments.

Graph Methods

DelaunayTriangulation.add_neighbour!Method
add_neighbour!(tri::Triangulation, u, v...)

Calls add_neighbour!(get_graph(tri), u, v...), adding each v into the neighbourhood of each u.

DelaunayTriangulation.each_vertexMethod
each_vertex(tri::Triangulation)

Returns each_vertex(get_graph(tri)), an iterator over all vertices present in the triangulation.

Warning

This iterator will include ghost vertices. If you want to exclude these, see each_solid_vertex. Alternatively, if you only want ghost vertices, see each_ghost_vertex.

Points Methods

DelaunayTriangulation.num_pointsMethod
num_points(tri::Triangulation)

Returns the number of points of tri.

Note

Note that this is just the size of get_points(tri), but if there are some missing points then this will not match the number in the triangulation itself. Use num_vertices to get the number of vertices in the triangulation.

DelaunayTriangulation.push_point!Method
push_point!(tri::Triangulation, x, y) = push_point!(get_points(tri), x, y)
push_point!(tri::Triangulation, p) = push_point!(get_points(tri), p)

Pushes the point p = (x, y) into get_points(tri0.

Note

This does not add the point p into the triangulation itself. See add_point! for this.

DelaunayTriangulation.each_solid_vertexFunction
each_solid_vertex(tri::Triangulation)

Returns an iterator over the solid vertices in the triangulation tri, i.e. those that are not ghost vertices.

DelaunayTriangulation.num_solid_verticesFunction
num_solid_vertices(tri::Triangulation)

Returns the number of solid vertices of tri.

num_solid_vertices(stats::TriangulationStatistics)

Returns the numsolidvertices field from the triangulation statistics stats.

DelaunayTriangulation.num_ghost_verticesFunction
num_ghost_vertices(tri::Trianngulation)

Returns the number of ghost vertices of tri.

num_ghost_vertices(stats::TriangulationStatistics)

Returns the numghostvertices field from the triangulation statistics stats.

Triangles Methods

DelaunayTriangulation.triangle_typeMethod
triangle_type(tri::Triangulation)

Returns the type used for representing triangles in tri.

triangle_type(vor::VoronoiTessellation)

Returns the type of the triangles used in the triangulation underlying the Voronoi tessellation vor.

DelaunayTriangulation.num_trianglesMethod
num_triangles(tri::Triangulation) = num_triangles(get_triangles(tri))

Returns the number of triangles in tri, including ghost triangles.

Predicates Methods

DelaunayTriangulation.is_boundary_edgeMethod
is_boundary_edge(tri::Triangulation, ij)
is_boundary_edge(tri::Triangulation, i, j)

Returns is_boundary_edge(ij, get_adjacent(tri)) or is_boundary_edge(i, j, get_adjacent(tri)), respectively, testing if ij = (i, j) belongs to the boundary of the triangulation tri, i.e. (i, j) adjoins a ghost vertex.

DelaunayTriangulation.is_boundary_triangleMethod
is_boundary_triangle(tri::Triangulation, i, j, k)
is_boundary_triangle(tri::Triangulation, T)

Returns is_boundary_triangle(i, j, k, get_adjacent(tri)) or is_boundary_triangle(T, get_adjacent(tri)), testing if at least one of the edges (u, v) of T = (i, j, k) satisfies is_boundary_edge(tri, u, v).

DelaunayTriangulation.point_position_relative_to_circumcircleMethod
point_position_relative_to_circumcircle(tri::Triangulation, i, j, k, ℓ)
point_position_relative_to_circumcircle(tri::Triangulation, T, ℓ)

Computes the position of the th point of tri relative to the circumcircle of the triangle T = (i, j, k) in tri.

DelaunayTriangulation.point_position_relative_to_triangleMethod
point_position_relative_to_triangle(tri::Triangulation, i, j, k, u)
point_position_relative_to_triangle(tri::Triangulation, T, u)

Computes the position of the uth point of tri relative to the triangle T = (i, j, k) in tri.

DelaunayTriangulation.triangle_line_segment_intersectionMethod
triangle_line_segment_intersection(tri::Triangulation, i, j, k, u, v)
triangle_line_segment_intersection(tri::Triangulation, T, e)

Computes the type of the intersection of the triangle T = (i, j, k) in tri with the line segment e = (u, v).

DelaunayTriangulation.is_outer_boundary_indexMethod
is_outer_boundary_index(tri::Triangulation, i)

Returns true if is_boundary_index(i) and this boundary index i corresponds to the outermost boundary of tri. Returns false otherwise.

DelaunayTriangulation.edge_existsMethod
edge_exists(tri::Triangulation, i, j)
edge_exists(tri::Triangulation, ij)

Returns true if the edge (i, j) exists in tri, and false otherwise.

DelaunayTriangulation.is_legalMethod
is_legal(tri::Triangulation, i, j)

Returns true if the edge (i, j) is legal in tri, and false otherwise. We also define constrained edges to be legal, as are ghost edges.

DelaunayTriangulation.find_edgeMethod
find_edge(tri::Triangulation, T, ℓ)

Given a point that is on an edge of the triangle T in tri, returns the edge that is on.

Representative Points Methods

DelaunayTriangulation.compute_representative_points!Method
compute_representative_points!(tri::Triangulation; use_convex_hull=!has_multiple_segments(tri) && num_boundary_edges(get_boundary_nodes(tri)) == 0)

Updates get_representative_point_list(tri) to match the current position of the boundaries. If there are no boundary nodes, use_convex_hull instead represents them using the indices of the convex hull.

DelaunayTriangulation.update_centroid_after_addition!Method
update_centroid_after_addition!(tri::Triangulation, i, p)

After the point p has been added into tri, this updates the ith representative point in get_representative_point_list(tri), treating it as if it were a centroid.

DelaunayTriangulation.update_centroid_after_deletion!Method
update_centroid_after_deletion!(tri::Triangulation, i, p)

After the point p has been deleted from tri, this updates the ith representative point in get_representative_point_list(tri), treating it as if it were a centroid.