Utility Functions
This section lists some of the internal utility functions, or other miscellaneous functions, used in this package.
DelaunayTriangulation.number_type — Function
number_type(x) -> DataTypeGiven a container x, returns the number type used for storing coordinates.
Examples
julia> using DelaunayTriangulation
julia> DelaunayTriangulation.number_type([1, 2, 3])
Int64
julia> DelaunayTriangulation.number_type((1, 2, 3))
Int64
julia> DelaunayTriangulation.number_type([1.0 2.0 3.0; 4.0 5.0 6.0])
Float64
julia> DelaunayTriangulation.number_type([[[1, 2, 3, 4, 5, 1]], [[6, 8, 9], [9, 10, 11], [11, 12, 6]]])
Int64
julia> DelaunayTriangulation.number_type((1.0f0, 2.0f0))
Float32
julia> DelaunayTriangulation.number_type(Vector{Float64})
Float64
julia> DelaunayTriangulation.number_type(Vector{Vector{Float64}})
Float64
julia> DelaunayTriangulation.number_type(NTuple{2, Float64})
Float64DelaunayTriangulation._to_val — Method
_to_val(v) -> ValWraps v in a Val, or if v isa Val simply returns v.
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.
DelaunayTriangulation.check_absolute_precision — Method
check_absolute_precision(x, y) -> BoolReturns true if abs(x - y) is less than or equal to sqrt(eps(Float64)).
DelaunayTriangulation.check_precision — Method
check_precision(x) -> BoolReturns true if abs(x) is less than or equal to sqrt(eps(number_type(eps))).
DelaunayTriangulation.check_ratio_precision — Method
check_ratio_precision(x, y) -> BoolReturns true if abs(x/y) is bounded between 0.99 and 1.01.
DelaunayTriangulation.check_relative_precision — Method
check_relative_precision(x, y) -> BoolReturns true if abs(x - y)/max(abs(x), abs(y)) is less than or equal to sqrt(eps(Float64)).
DelaunayTriangulation.choose_uvw — Method
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).
DelaunayTriangulation.circular_equality — Method
circular_equality(A, B, by=isequal) -> BoolCompares the two circular vectors A and B for equality up to circular rotation, using by to compare individual elements.
DelaunayTriangulation.dist — Method
dist(p, q) -> NumberAssuming p and q are two-dimensional, computes the Euclidean distance between p and q.
DelaunayTriangulation.dist — Method
dist(tri::Triangulation, p) -> NumberGiven 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, thenpis inside the triangulation.δ < 0: If the returned distance is negative, thenpis outside the triangulation.δ = 0: If the returned distance is zero, thenpis on the boundary of the triangulation.
Where we say distance, we are referring to the distance from p to the boundary of the triangulation.
DelaunayTriangulation.dist_sqr — Method
dist_sqr(p, q) -> NumberAssuming p and q are two-dimensional, computes the square of the Euclidean distance between p and q.
DelaunayTriangulation.edge_length — Function
edge_length(tri::Triangulation, u, v) -> Number
edge_length(tri::Triangulation, e) -> NumberComputes the length of the edge e = (u, v).
DelaunayTriangulation.edge_length_sqr — Method
edge_length_sqr(tri::Triangulation, u, v) -> Number
edge_length_sqr(tri::Triangulation, e) -> NumberComputes the square of the length of the edge e = (u, v).
DelaunayTriangulation.eval_fnc_at_het_tuple_element — Method
eval_fnc_at_het_tuple_element(f, tup, idx)Evaluates f(tup[idx]) in a type-stable way. If idx > length(tup), then f is evaluated on the last element of tup. If length(tup) > 32, then the function is not type-stable; note that, in this case, idx > length(tup) leads to a BoundsError.
DelaunayTriangulation.eval_fnc_at_het_tuple_element_with_arg — Method
eval_fnc_at_het_tuple_element_with_arg(f, tup, arg, idx)Evaluates f(tup[idx], arg...) in a type-stable way. If idx > length(tup), then f is evaluated on the last element of tup. If length(tup) > 32, then the function is not type-stable; note that, in this case, idx > length(tup) leads to a BoundsError.
DelaunayTriangulation.eval_fnc_at_het_tuple_element_with_arg_and_prearg — Method
eval_fnc_at_het_tuple_element_with_arg_and_prearg(f, tup, prearg, arg, idx)Evaluates f(prearg, tup[idx], arg...) in a type-stable way. If idx > length(tup), then f is evaluated on the last element of tup. If length(tup) > 32, then the function is not type-stable; note that, in this case, idx > length(tup) leads to a BoundsError.
DelaunayTriangulation.eval_fnc_at_het_tuple_two_elements — Method
eval_fnc_at_het_tuple_two_elements(f, tup, idx1, idx2)Evaluates f(tup[idx1], tup[idx2]) in a type-stable way. If length(tup) > 32, then the function is not type-stable.
DelaunayTriangulation.eval_fnc_in_het_tuple — Method
eval_fnc_in_het_tuple(tup, arg, idx)Evaluates tup[idx](arg...) in a type-stable way. If idx > length(tup), then tup[end](arg...) is evaluated.
DelaunayTriangulation.get_area — Method
get_area(tri::Triangulation) -> NumberReturns the area of tri.
DelaunayTriangulation.get_boundary_chain — Method
get_boundary_chain(tri::Triangulation, i, j) -> EdgesGiven 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.
DelaunayTriangulation.get_ghost_vertex — Function
get_ghost_vertex(i, j, k) -> Vertex
get_ghost_vertex(i, j) -> VertexGiven 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)
-5DelaunayTriangulation.get_ordinal_suffix — Method
get_ordinal_suffix(i) -> StringReturns 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"DelaunayTriangulation.is_circular — Method
is_circular(A) -> BoolTests if A is circular, meaning that A[begin] == A[end].
DelaunayTriangulation.is_true — Function
is_true(b) -> BoolReturns 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))
falseDelaunayTriangulation.iterated_neighbourhood — Method
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$.
DelaunayTriangulation.midpoint — Method
midpoint(p, q) -> Number or NTuple{2, Number}Assuming p and q are either both numbers are both 2-vectors, computes their average.
DelaunayTriangulation.midpoint — Method
midpoint(tri::Triangulation, u, v) -> NTuple{2, Number}
midpoint(tri::Triangulation, e) -> NTuple{2, Number}Computes the midpoint of e = (u, v).
DelaunayTriangulation.nextindex_circular — Method
nextindex_circular(C, i) -> IntegerReturns the next index after i in the circular vector C.
DelaunayTriangulation.norm — Method
norm(p) -> NumberAssuming p is two-dimensional, computes the Euclidean norm of p.
DelaunayTriangulation.norm_sqr — Method
norm_sqr(p) -> NumberAssuming p is two-dimensional, computes the square of the Euclidean norm of p.
DelaunayTriangulation.previndex_circular — Method
previndex_circular(C, i) -> IntegerReturns the previous index before i in the circular vector C.
DelaunayTriangulation.replace_boundary_triangle_with_ghost_triangle — Method
replace_boundary_triangle_with_ghost_triangle(tri::Triangulation, V) -> TriangleGiven a boundary triangle V of tri, returns the adjacent ghost triangle. Note that for triangles in a corner of a domain, like a lattice triangulation, there are two choices of ghost triangle.
DelaunayTriangulation.replace_ghost_triangle_with_boundary_triangle — Method
replace_ghost_triangle_with_boundary_triangle(tri::Triangulation, V) -> TriangleGiven a ghost triangle V of tri, returns the adjacent boundary triangle.
DelaunayTriangulation.self_eval — Method
self_eval(f, args...)Evaluates f(args...).
DelaunayTriangulation.uniquetol — Method
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.
DelaunayTriangulation.angle_between — Method
angle_between(p, q) -> NumberReturns the angle between the vectors p and q in radians, treating q as the base. See this article. The returned angle is in [0, 2π).
DelaunayTriangulation.classify_and_compute_segment_intersection — Method
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: ACertificateindicating the intersection type.cert_c: ACertificateindicating the position ofcrelative to the line through(a, b).cert_d: ACertificateindicating the position ofdrelative to the line through(a, b).p: The intersection point ifcertisCert.SingleorCert.Touching, and(NaN, NaN)otherwise.
DelaunayTriangulation.distance_to_polygon — Method
distance_to_polygon(q, points, boundary_nodes) -> NumberGiven 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.
DelaunayTriangulation.get_distance_to_plane — Method
get_distance_to_plane(a, b, c, p) -> NumberReturns the distance from the point p to the plane defined by the points (a, b, c). The distance is positive if p is above the plane.
DelaunayTriangulation.get_plane_through_three_points — Method
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*}\]
DelaunayTriangulation.get_steepest_descent_direction — Method
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 + δ = 0is the plane, then the steepest descent direction is (α, β)/γ. The returned value is given by (x, y) = sign(γ)(α, β).
See also get_plane_through_three_points.
DelaunayTriangulation.get_vertical_distance_to_plane — Method
get_vertical_distance_to_plane(a, b, c, p) -> NumberReturns the vertical distance from the point p to the plane defined by the points (a, b, c). The distance is positive if p is above the plane.
DelaunayTriangulation.intersection_of_edge_and_bisector_ray — Method
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: ACertificateindicating the position ofcrelative to the line through(a, b).p: The intersection point ifcintersects 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.
DelaunayTriangulation.intersection_of_ray_with_bounding_box — Method
intersection_of_ray_with_bounding_box(p, q, a, b, c, d) -> NTuple{2, Number}Compute the intersection of the ray emanating from p and passing through q with the box [a, b] × [c, d]. It is assumed that p is inside of the box.
DelaunayTriangulation.polygon_bounds — Function
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.
DelaunayTriangulation.polygon_features — Method
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.
DelaunayTriangulation.project_onto_line — Method
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].
DelaunayTriangulation.segment_intersection_coordinates — Method
segment_intersection_coordinates(a, b, c, d) -> NTuple{2, Number}Given two segments (a, b) and (c, d) that are assumed to intersect, computes the coordinates of the intersection.
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].
DelaunayTriangulation.squared_distance_to_segment — Method
squared_distance_to_segment(x₁, y₁, x₂, y₂, x, y) -> NumberGiven a line segment (x₁, y₁) → (x₂, y₂) and a query point (x, y), returns the squared distance from (x, y) to the line segment.
DelaunayTriangulation.angle_is_acute_predicate — Method
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:∠prqis acute.0:∠prqis a right angle.-1:∠prqis 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.
DelaunayTriangulation.det — Method
det(a, b, c, d)Computes ExactPredicates.det(a, b, c, d), i.e. returns a*d - b*c.
DelaunayTriangulation.det — Method
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}\]
DelaunayTriangulation.ext — Method
ext(ux, uy, vx, vy)Computes ExactPredicates.ext((ux, uy), (vx, vy)), i.e. returns ux * vy - uy * vx.
DelaunayTriangulation.incircle_predicate — Method
incircle_predicate([kernel::AbstractPredicateKernel,] a, b, c, p; cache = nothing) -> IntegerAssuming that (a, b, c) is a positively oriented triangle, returns
1: Ifpis inside the circle defined by(a, b, c).0: Ifpis on the circle defined by(a, b, c).-1: Ifpis 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.
DelaunayTriangulation.inp — Method
inp(ux, uy, vx, vy)Computes ExactPredicates.inp((ux, uy), (vx, vy)), i.e. returns ux * vx + uy * vy.
DelaunayTriangulation.meet_predicate — Method
meet_predicate([kernel::AbstractPredicateKernel], p, q, a, b) -> IntegerReturns
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.
DelaunayTriangulation.opposite_signs — Method
opposite_signs(x, y) -> BoolFrom ExactPredicates.jl, returns true if x and y have opposite signs, and false otherwise. Assumes that x and y are in [-1, 0, 1].
DelaunayTriangulation.orient_predicate — Method
orient_predicate([kernel::AbstractPredicateKernel,] p, q, r) -> IntegerReturns
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.
DelaunayTriangulation.orient_predicate — Method
orient_predicate([kernel::AbstractPredicateKernel,] p, q, r; cache = nothing) -> IntegerReturns
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.
DelaunayTriangulation.parallelorder_predicate — Method
parallelorder_predicate([kernel::AbstractPredicateKernel,] a, b, p, q) -> IntegerReturns
1:qis closer to the line(a, b)thanp.0:pandqare equidistant from the line(a, b).-1:pis closer to the line(a, b)thanq.
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.
DelaunayTriangulation.sameside_predicate — Method
sameside_predicate(a, b, p) -> IntegerAssuming all of a, b, p are collinear, returns
1:aandbare on the same side ofpon the line.0:a == porb == p.-1:aandbare on different sides ofpon the line.
DelaunayTriangulation.sgn — Method
sgn(x) -> IntReturns Int(sign(x)).
DelaunayTriangulation.is_boundary_edge — Method
is_boundary_edge(tri::Triangulation, ij) -> Bool
is_boundary_edge(tri::Triangulation, i, j) -> BoolTests if the edge (i, j) is a boundary edge of tri, meaning (j, i) adjoins a ghost vertex.
DelaunayTriangulation.is_boundary_node — Method
is_boundary_node(tri::Triangulation, i) -> (Bool, Vertex)Tests if the vertex i is a boundary node of tri.
Arguments
tri::Triangulation: TheTriangulation.i: The vertex to test.
Outputs
flag:trueifiis a boundary node, andfalseotherwise.g: Either the ghost vertex corresponding with the section thatilives on ifflagis true, or 0 otherwise.
DelaunayTriangulation.is_boundary_triangle — Method
is_boundary_triangle(tri::Triangulation, T) -> Bool
is_boundary_triangle(tri::Triangulation, i, j, k) -> BoolReturns true if the triangle T = (i, j, k) of tri has an edge on the boundary, and false otherwise.
DelaunayTriangulation.is_exterior_boundary_node — Method
is_exterior_boundary_node(tri::Triangulation, i) -> BoolTests if the vertex i is an exterior boundary node of tri.
DelaunayTriangulation.is_exterior_ghost_edge — Method
is_exterior_ghost_edge(tri::Triangulation, i, j) -> BoolTests if the edge (i, j) is an exterior ghost edge of tri.
See also is_exterior_ghost_vertex.
DelaunayTriangulation.is_exterior_ghost_triangle — Method
is_exterior_ghost_triangle(tri::Triangulation, i, j, k) -> BoolTests if the triangle (i, j, k) is an exterior ghost triangle of tri.
See also is_exterior_ghost_vertex.
DelaunayTriangulation.is_ghost_edge — Method
is_ghost_edge(ij) -> Bool
is_ghost_edge(i, j) -> BoolTests if the edge (i, j) is a ghost edge, meaning i or j is a ghost vertex.
DelaunayTriangulation.is_ghost_triangle — Method
is_ghost_triangle(T) -> Bool
is_ghost_triangle(i, j, k) -> BoolTests if T = (i, j, k) is a ghost triangle, meaning i, j or k is a ghost vertex.
DelaunayTriangulation.is_ghost_vertex — Method
is_ghost_vertex(i) -> BoolTests if i is a ghost vertex, meaning i ≤ -1.
DelaunayTriangulation.convert_certificate — Function
convert_certificate(cert::I, Cert1, Cert2, Cert3) -> CertificateGiven cert ∈ (-1, 0, 1), return Cert1, Cert2 or Cert3 depending on if cert == -1, cert == 0 or cert == 1, respectively.
DelaunayTriangulation.DefaultAdjacentValue — Constant
DefaultAdjacentValue = 0Default value used for representing an empty result from an adjacency query.
DelaunayTriangulation.𝒢 — Constant
𝒢 = GhostVertexAlias for GhostVertex.
DelaunayTriangulation.GhostVertex — Constant
GhostVertex = -1Number 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/.
DelaunayTriangulation.ε — Function
ε(x) = sqrt(eps(number_type(x)))Number used as a tolerance in certain functions, e.g. for mesh refinement when using check_precision to avoid degenerate circumcenters.
DelaunayTriangulation.∅ — Constant
∅ = DefaultAdjacentValueAlias for DefaultAdjacentValue.
DelaunayTriangulation.fix_orient3_cache — Function
fix_orient3_cache(tri::Triangulation, cache)Returns cache if validate_orient3_cache(tri, cache) is true, otherwise returns nothing.
DelaunayTriangulation.fix_incircle_cache — Function
fix_incircle_cache(tri::Triangulation, cache)Returns cache if validate_incircle_cache(tri, cache) is true, otherwise returns nothing.
DelaunayTriangulation.validate_orient3_cache — Function
validate_orient3_cache(tri::Triangulation, cache) -> BoolChecks 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.
DelaunayTriangulation.validate_incircle_cache — Function
validate_incircle_cache(tri::Triangulation, cache) -> BoolChecks 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.