Utility Functions

This section lists some of the internal utility functions, or other miscellaneous functions, used in this package.

DelaunayTriangulation.adjust_θMethod
adjust_θ(θ₁, θ₂, positive) -> (Number, Number)

Given two angles θ₁ and θ₂ in radians, adjusts the angles to new angles θ₁′, θ₂′ so that θ₁′ ≤ θ₂′ if positive is true, and θ₁′ ≥ θ₂′ if positive is false.

source
DelaunayTriangulation.choose_uvwMethod
choose_uvw(e1, e2, e3, u, v, w) -> (Vertex, Vertex, Vertex)

Choose values for (u, v, w) based on the Booleans (e1, e2, e3), assuming only one is true. The three cases are:

  • If e1, returns (u, v, w).
  • If e2, returns (v, w, u).
  • If e3, returns (w, u, v).
source
DelaunayTriangulation.distMethod
dist(tri::Triangulation, p) -> Number

Given a point p, returns the distance from p to the triangulation, using the conventions from distance_to_polygon:

  • δ > 0: If the returned distance is positive, then p is inside the triangulation.
  • δ < 0: If the returned distance is negative, then p is outside the triangulation.
  • δ = 0: If the returned distance is zero, then p is on the boundary of the triangulation.

Where we say distance, we are referring to the distance from p to the boundary of the triangulation.

source
DelaunayTriangulation.get_boundary_chainMethod
get_boundary_chain(tri::Triangulation, i, j) -> Edges

Given two boundary vertices i and j on a boundary with ghost vertex ghost_vertex, walks counter-clockwise from i to j along the boundary and returns the collection of all vertices encountered in counter-clockwise order.

source
DelaunayTriangulation.get_ghost_vertexFunction
get_ghost_vertex(i, j, k) -> Vertex
get_ghost_vertex(i, j) -> Vertex

Given three vertices i, j, and k, returns the ghost vertex among them. If none of them are ghost vertices, returns k. The two-argument version is equivalent to get_ghost_vertex(i, j, j).

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.get_ghost_vertex(1, 7, -2)
-2

julia> DelaunayTriangulation.get_ghost_vertex(-1, 2, 3)
-1

julia> DelaunayTriangulation.get_ghost_vertex(1, 5, 10)
10

julia> DelaunayTriangulation.get_ghost_vertex(1, -1)
-1

julia> DelaunayTriangulation.get_ghost_vertex(-5, 2)
-5
source
DelaunayTriangulation.get_ordinal_suffixMethod
get_ordinal_suffix(i) -> String

Returns the ordinal suffix for the integer i.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.get_ordinal_suffix(1)
"st"

julia> DelaunayTriangulation.get_ordinal_suffix(2)
"nd"

julia> DelaunayTriangulation.get_ordinal_suffix(3)
"rd"

julia> DelaunayTriangulation.get_ordinal_suffix(4)
"th"

julia> DelaunayTriangulation.get_ordinal_suffix(5)
"th"

julia> DelaunayTriangulation.get_ordinal_suffix(6)
"th"

julia> DelaunayTriangulation.get_ordinal_suffix(11)
"th"

julia> DelaunayTriangulation.get_ordinal_suffix(15)
"th"

julia> DelaunayTriangulation.get_ordinal_suffix(100)
"th"
source
DelaunayTriangulation.is_trueFunction
is_true(b) -> Bool

Returns b represents a true value, and false otherwise.

Examples

julia> using DelaunayTriangulation

julia> DelaunayTriangulation.is_true(true)
true

julia> DelaunayTriangulation.is_true(false)
false

julia> DelaunayTriangulation.is_true(Val(true))
true

julia> DelaunayTriangulation.is_true(Val(false))
false
source
DelaunayTriangulation.iterated_neighbourhoodMethod
iterated_neighbourhood(tri::Triangulation, i, d) -> Set{Vertex}

Returns the set of vertices of tri in the iterated neighbourhood of the vertex i of depth d.

Extended help

The $d$-times iterated neighbourhood is defined by

\[N_i^d = \bigcup_{j \in N_i^{d-1}} N_j \setminus \{i\},\]

where $N_i^1 = N_i$ is the set of neighbours of $i$.

source
DelaunayTriangulation.midpointMethod
midpoint(tri::Triangulation, u, v) -> NTuple{2, Number}
midpoint(tri::Triangulation, e) -> NTuple{2, Number}

Computes the midpoint of e = (u, v).

source
DelaunayTriangulation.uniquetolMethod
uniquetol(A::Vector{Float64}; tol=1e-12) -> Vector{Float64}

Returns the unique elements of A up to tolerance tol. We say that two values x and y are within tolerance if abs(u - v) ≤ M*tol, where M = maximum(abs.(A)). It is assumed that A is sorted - this is NOT checked.

source
DelaunayTriangulation.classify_and_compute_segment_intersectionMethod
classify_and_compute_segment_intersection([kernel::AbstractPredicateKernel,] a, b, c, d) -> (Certificate, Certificate, Certificate, NTuple{2, Number})

Given two line segments (a, b) and (c, d), classifies the intersection of the two segments. The returned value is (cert, cert_c, cert_d, p), where:

  • cert: A Certificate indicating the intersection type.
  • cert_c: A Certificate indicating the position of c relative to the line through (a, b).
  • cert_d: A Certificate indicating the position of d relative to the line through (a, b).
  • p: The intersection point if cert is Cert.Single or Cert.Touching, and (NaN, NaN) otherwise.
source
DelaunayTriangulation.distance_to_polygonMethod
distance_to_polygon(q, points, boundary_nodes) -> Number

Given a query point q and a polygon defined by (points, boundary_nodes), returns the signed distance from q to the polygon. The boundary_nodes must match the specification in the documentation and in check_args.

See also dist.

source
DelaunayTriangulation.get_plane_through_three_pointsMethod
get_plane_through_three_points(a, b, c) -> NTuple{4, Number}

Given three points (a, b, c) in ℝ³ represented as Tuples, computes the equation of the plane through the points. The result is given in the form (α, β, γ, δ), so that the plane is given by

αx + βy + γz + δ = 0.

Extended help

The equation of the plane is computed by expanding the equation

\[\det \begin{bmatrix} x & y & z & 1 \\ a_x & a_y & a_z & 1 \\ b_x & b_y & b_z & 1 \\ c_x & c_y & c_z & 1 \end{bmatrix} = 0.\]

From this, we find:

\[\begin{align*} \alpha &= a_y b_z - a_z b_y - a_y c_z + a_z c_y + b_y c_z - b_z c_y, \\ \beta &= a_z b_x - a_x b_z + a_x c_z - a_z c_x - b_x c_z + b_z c_x, \\ \gamma &= a_x b_y - a_y b_x - a_x c_y + a_y c_x + b_x c_y - b_y c_x, \\ \delta &= a_x b_z c_y - a_x b_y c_z + a_y b_x c_z - a_y b_z c_x - a_z b_x c_y + a_z b_y c_x. \end{align*}\]

source
DelaunayTriangulation.get_steepest_descent_directionMethod
get_steepest_descent_direction(a, b, c) -> NTuple{2, Number}

Given three points in ℝ³ defining a plane, returns the direction (x, y) of the steepest descent along the plane. In particular, if

αx + βy + γz + δ = 0

is the plane, then the steepest descent direction is (α, β)/γ. The returned value is given by (x, y) = sign(γ)(α, β).

See also get_plane_through_three_points.

source
DelaunayTriangulation.intersection_of_edge_and_bisector_rayMethod
intersection_of_edge_and_bisector_ray([kernel::AbstractPredicateKernel=AdaptiveKernel(),] a, b, c; project=false) -> (Certificate, NTuple{2, Number})

Given an edge (a, b) and a ray emanating from c perpendicular with the edge and collinear with its midpoint (or, if project=true, the projection of c onto the edge), tests if c intersects the edge. The returned value is (cert, p), where:

  • cert: A Certificate indicating the position of c relative to the line through (a, b).
  • p: The intersection point if c intersects the edge, (NaN, NaN) otherwise.

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

source
DelaunayTriangulation.polygon_boundsFunction
polygon_bounds(points, boundary_nodes, check_all_curves=Val(false)) -> (Number, Number, Number, Number)

Computes the bounding box of the polygon defined by (points, boundary_nodes). The boundary_nodes must match the specification in the documentation and in check_args. If check_all_curves is true, then the bounding box of the union of all curves of the polygon is computed instead of just the first curve.

source
DelaunayTriangulation.polygon_featuresMethod
polygon_features(points, boundary_nodes) -> (Number, NTuple{2, Number})

Computes the signed area and centroid of the polygon defined by (points, boundary_nodes). The boundary_nodes must match the specification in the documentation and in check_args.

source
DelaunayTriangulation.project_onto_lineMethod
project_onto_line(p, q, r) -> NTuple{2, Number}

Projects the point r onto the line through p and q. It is possible that the projected point is not on the line segment [p, q].

source
DelaunayTriangulation.sort_convex_polygon!Method
sort_convex_polygon!(vertices, points)

Sorts the vertices of a convex polygon in counter-clockwise order. The polygon is defined by (points, vertices), and the vertices are sorted in-place. It is assumed that the vertices are not circular, i.e. vertices[begin] ≠ vertices[end].

source
DelaunayTriangulation.angle_is_acute_predicateMethod
angle_is_acute([kernel::AbstractPredicateKernel,] p, q, r)

Tests if the angle opposite (p, q) in the triangle (p, q, r), meaning ∠prq, is acute, returning:

  • 1: ∠prq is acute.
  • 0: ∠prq is a right angle.
  • -1: ∠prq is obtuse.

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

source
DelaunayTriangulation.detMethod
det(a, b, c, d, e, f, g, h, i)

Computes ExactPredicates.det(a, b, c, d, e, f, g, h, i), i.e. returns the determinant of

\[\det \begin{bmatrix} a & b & c \\ d & e & f \\ g & h & i \end{bmatrix}\]

source
DelaunayTriangulation.incircle_predicateMethod
incircle_predicate([kernel::AbstractPredicateKernel,] a, b, c, p; cache = nothing) -> Integer

Assuming that (a, b, c) is a positively oriented triangle, returns

  • 1: If p is inside the circle defined by (a, b, c).
  • 0: If p is on the circle defined by (a, b, c).
  • -1: If p is outside the circle defined by (a, b, c).

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

The optional cache keyword argument can be used for preallocating memory for intermediate results, passing the argument from AdaptivePredicates.incircleadapt_cache(T), where T is the number type of the input points. If nothing is passed, no cache is used. This is only needed if an AdaptiveKernel() is used.

source
DelaunayTriangulation.meet_predicateMethod
meet_predicate([kernel::AbstractPredicateKernel], p, q, a, b) -> Integer

Returns

  • 1: The open line segments (p, q) and (a, b) meet in a single point.
  • 0: The closed line segments [p, q] and [a, b] meet in one or several points.
  • -1: Otherwise.

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

source
DelaunayTriangulation.opposite_signsMethod
opposite_signs(x, y) -> Bool

From ExactPredicates.jl, returns true if x and y have opposite signs, and false otherwise. Assumes that x and y are in [-1, 0, 1].

source
DelaunayTriangulation.orient_predicateMethod
orient_predicate([kernel::AbstractPredicateKernel,] p, q, r) -> Integer

Returns

  • 1: (p, q, r) is positively oriented.
  • 0: (p, q, r) is collinear / degenerate.
  • -1: (p, q, r) is negatively oriented.

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

source
DelaunayTriangulation.orient_predicateMethod
orient_predicate([kernel::AbstractPredicateKernel,] p, q, r; cache = nothing) -> Integer

Returns

  • 1: (p, q, r, s) is positively oriented.
  • 0: (p, q, r, s) is collinear / degenerate.
  • -1: (p, q, r, s) is negatively oriented.

Here, a positively oriented tetrahedron (p, q, r, s) takes the form

                               z.
                             .
                           ,/
                         s
                       ,/|'\
                     ,/  |  '\
                   ,/    '.   '\
                 ,/       |     '\                 
               ,/         |       '\              
              p-----------'.--------q --> x
               '\.         |      ,/              
                  '\.      |    ,/                 
                     '\.   '. ,/    
                        '\. |/      
                           'r       
                             '\.

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

The optional cache keyword argument can be used for preallocating memory for intermediate results, passing the argument from AdaptivePredicates.orient3adapt_cache(T), where T is the number type of the input points. If nothing is passed, no cache is used. This is only needed if an AdaptiveKernel() is used.

source
DelaunayTriangulation.parallelorder_predicateMethod
parallelorder_predicate([kernel::AbstractPredicateKernel,] a, b, p, q) -> Integer

Returns

  • 1: q is closer to the line (a, b) than p.
  • 0: p and q are equidistant from the line (a, b).
  • -1: p is closer to the line (a, b) than q.

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

source
DelaunayTriangulation.sameside_predicateMethod
sameside_predicate(a, b, p) -> Integer

Assuming all of a, b, p are collinear, returns

  • 1: a and b are on the same side of p on the line.
  • 0: a == p or b == p.
  • -1: a and b are on different sides of p on the line.
Note

The difference in the argument order to ExactPredicates.jl is to match the convention that the main point being tested is the last argument.

source
DelaunayTriangulation.is_boundary_edgeMethod
is_boundary_edge(tri::Triangulation, ij) -> Bool
is_boundary_edge(tri::Triangulation, i, j) -> Bool

Tests if the edge (i, j) is a boundary edge of tri, meaning (j, i) adjoins a ghost vertex.

source
DelaunayTriangulation.is_boundary_nodeMethod
is_boundary_node(tri::Triangulation, i) -> (Bool, Vertex)

Tests if the vertex i is a boundary node of tri.

Arguments

Outputs

  • flag: true if i is a boundary node, and false otherwise.
  • g: Either the ghost vertex corresponding with the section that i lives on if flag is true, or 0 otherwise.
source
DelaunayTriangulation.is_boundary_triangleMethod
is_boundary_triangle(tri::Triangulation, T) -> Bool
is_boundary_triangle(tri::Triangulation, i, j, k) -> Bool

Returns true if the triangle T = (i, j, k) of tri has an edge on the boundary, and false otherwise.

source
DelaunayTriangulation.convert_certificateFunction
convert_certificate(cert::I, Cert1, Cert2, Cert3) -> Certificate

Given cert ∈ (-1, 0, 1), return Cert1, Cert2 or Cert3 depending on if cert == -1, cert == 0 or cert == 1, respectively.

source
DelaunayTriangulation.GhostVertexConstant
GhostVertex = -1

Number used for representing initial ghost vertices. All other ghost vertices are derived from subtracting from this number. See https://juliageometry.github.io/DelaunayTriangulation.jl/stable/manual/ghost_triangles/.

source
DelaunayTriangulation.validate_orient3_cacheFunction
validate_orient3_cache(tri::Triangulation, cache) -> Bool

Checks if the cache cache is of the correct type for cmoputing orient3 predicates for the triangulation tri. If isnothing(cache) or the cache is of the correct type, then true is returned. Otherwise, false is returned.

source
DelaunayTriangulation.validate_incircle_cacheFunction
validate_incircle_cache(tri::Triangulation, cache) -> Bool

Checks if the cache cache is of the correct type for computing incircle predicates for the triangulation tri. If isnothing(cache) or the cache is of the correct type, then true is returned. Otherwise, false is returned.

source