Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix issue with clipped Voronoi tessellation for a single right-angled triangle #207

Merged
merged 4 commits into from
Oct 19, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 3 additions & 0 deletions NEWS.md
Original file line number Diff line number Diff line change
@@ -1,5 +1,8 @@
# Changelog

## 1.6.1
- Fix issue with clipping Voronoi tessellation dual to a single right-angled triangle. See [#207](https://github.com/JuliaGeometry/DelaunayTriangulation.jl/pull/207)

## 1.6.0
- Define `reverse` for `AbstractParametricCurve`s, making it easier to reverse the orientation of a curve. See [#195](https://github.com/JuliaGeometry/DelaunayTriangulation.jl/pull/195).
- Fixed an issue with `LineSegment` not returning the exact endpoints at `t=1`, which can be problematic when joining boundary nodes. This has been fixed. See [#195](https://github.com/JuliaGeometry/DelaunayTriangulation.jl/pull/195).
Expand Down
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "DelaunayTriangulation"
uuid = "927a84f5-c5f4-47a5-9785-b46e178433df"
authors = ["Daniel VandenHeuvel <[email protected]>"]
version = "1.6.0"
version = "1.6.1"

[deps]
AdaptivePredicates = "35492f91-a3bd-45ad-95db-fcad7dcfedb7"
Expand Down
2 changes: 1 addition & 1 deletion src/algorithms/voronoi/clipped_coordinates.jl
Original file line number Diff line number Diff line change
Expand Up @@ -151,7 +151,7 @@ function _get_ray(vorn::VoronoiTessellation, i, ghost_vertex, predicates::Abstra
m = midpoint(p, q)
end
mx, my = getxy(m)
if r == m # It's possible for the circumcenter to lie on the edge and exactly at the midpoint (e.g. [(0.0,1.0),(-1.0,2.0),(-2.0,-1.0)]). In this case, just rotate
if check_precision(dist(r, m)) # It's possible for the circumcenter to lie on the edge and exactly at the midpoint (e.g. [(0.0,1.0),(-1.0,2.0),(-2.0,-1.0)]). In this case, just rotate
dx, dy = qx - mx, qy - my
if is_right(point_position_relative_to_line(predicates, p, q, r))
rotated_dx, rotated_dy = dy, -dx
Expand Down
36 changes: 27 additions & 9 deletions test/voronoi/voronoi.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1384,7 +1384,7 @@ end
@test validate_tessellation(vorn, check_convex=false)
end

@testset "Clipping to a generic convex polygon" begin
@testset "Clipping to a generic convex polygon" begin
a, b, c, d, e, f, g, h, i = (-5.0, 7.0),
(-8.0, 6.0), (-7.0, 4.0), (-4.0, 3.0), (-6.0, 6.0),
(-1.0, 6.0), (-4.0, 6.0), (-3.0, 5.0), (-2.0, 8.0)
Expand Down Expand Up @@ -1482,11 +1482,11 @@ end
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[7]...), pt7, ⪧)
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[8]...), pt8, ⪧)
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[9]...), pt9, ⪧)
@test validate_tessellation(vorn; check_area = false)
@test validate_tessellation(vorn; check_area=false)
@test isempty(vorn.unbounded_polygons)
@test vorn.boundary_polygons == Set((4, 6, 2, 9, 3, 1))

vorn = voronoi(tri, clip=true, clip_polygon=clip_polygon, rng=rng, smooth = true)
vorn = voronoi(tri, clip=true, clip_polygon=clip_polygon, rng=rng, smooth=true)
pt1 = [
(-5.1590315116169645, 8.028359015649276)
(-2.0329830757675342, 9.155154224395774)
Expand Down Expand Up @@ -1551,7 +1551,7 @@ end
(-0.15539885081025462, 7.000135534420414)
(2.0, 7.000229455337406)
]
@test validate_tessellation(vorn; check_area = false)
@test validate_tessellation(vorn; check_area=false)
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[1]...), pt1, ⪧)
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[2]...), pt2, ⪧)
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[3]...), pt3, ⪧)
Expand All @@ -1562,7 +1562,7 @@ end
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[9]...), pt9, ⪧)
for index in each_polygon_index(vorn)
c = DT.get_centroid(vorn, index)
@test c ⪧ get_generator(vorn, index) rtol=1e-2 atol=1e-4
@test c ⪧ get_generator(vorn, index) rtol = 1e-2 atol = 1e-4
end

J, K, L, M = (-2.0, 12.0), (-8.0, 4.0), (-4.0, 2.0), (-2.0, 10.0)
Expand Down Expand Up @@ -1711,7 +1711,7 @@ end
@test DT.circular_equality(DT.get_polygon_point(vorn, vorn.polygons[9]...), pt9, ⪧)
for index in each_polygon_index(vorn)
c = DT.get_centroid(vorn, index)
@test c ⪧ get_generator(vorn, index) rtol=1e-2 atol=1e-4
@test c ⪧ get_generator(vorn, index) rtol = 1e-2 atol = 1e-4
end
end

Expand All @@ -1722,7 +1722,7 @@ end
@test vorn == vorn
@test vorn == vorn2
g1 = vorn.generators[1]
vorn.generators[1] = vorn.generators[2]
vorn.generators[1] = vorn.generators[2]
@test vorn != vorn2
vorn.generators[1] = g1
@test vorn == vorn2
Expand Down Expand Up @@ -1761,8 +1761,8 @@ end
@test vorn == vorn2
end

@testset "copy/deepcopy" begin
vorn = voronoi(triangulate(rand(2, 100)), clip = true, smooth = true)
@testset "copy/deepcopy" begin
vorn = voronoi(triangulate(rand(2, 100)), clip=true, smooth=true)
vorn1 = deepcopy(vorn)
vorn2 = copy(vorn)
@test typeof(vorn1) == typeof(vorn) == typeof(vorn2)
Expand All @@ -1773,4 +1773,22 @@ end
@test !(getfield(vorn, f) === getfield(_vorn, f))
end
end
end

@testset "Issue #206" begin
Random.seed!(123)
pts = [(0.1, 0.2), (1.1, 0.2), (0.1, 1.7)]
tri = triangulate(pts)
vorn = voronoi(tri, clip=true, clip_polygon=(
[(0.0, 0.0), (1.0, 0.0), (1.0, 1.0), (0.0, 1.0)], [1, 2, 3, 4, 1]
))
@test validate_tessellation(vorn; check_area=false)
A = 0.0
for i in each_polygon_index(vorn)
A += get_area(vorn, i)
end
@test A ≈ 1.0
vorn = voronoi(tri)
poly = get_polygon_coordinates(vorn, 3, (0.0, 1.0, 0.0, 1.0))
@test poly ⪧ [(0.675, 1.0), (0.0, 1.0), (0.0, 0.95), (0.6, 0.95), (0.675, 1.0)]
end
Loading