diff --git a/Project.toml b/Project.toml index 1de293763..bbbcd70e8 100644 --- a/Project.toml +++ b/Project.toml @@ -5,6 +5,7 @@ version = "0.53.0" [deps] Bessels = "0e736298-9ec6-45e8-9647-e4fc86a2fe38" +BinaryTrees = "289e92be-c05a-437b-9e67-5b0c799891f8" CircularArrays = "7a955b69-7140-5f4e-a0ed-f168c5e2e749" Colorfy = "03fe91ce-8ec6-4610-8e8d-e7491ccca690" CoordRefSystems = "b46f11dc-f210-4604-bfba-323c1ec968cb" @@ -22,8 +23,15 @@ TiledIteration = "06e1c1a7-607b-532d-9fad-de7d9aa2abac" TransformsBase = "28dd2a49-a57a-4bfb-84ca-1a49db9b96b8" Unitful = "1986cc42-f94f-5a68-af5c-568840ba703d" +[weakdeps] +Makie = "ee78f7c6-11fb-53f2-987a-cfe4a2b5a57a" + +[extensions] +MeshesMakieExt = "Makie" + [compat] Bessels = "0.2" +BinaryTrees = "0.1.2" CircularArrays = "1.3" Colorfy = "1.1" CoordRefSystems = "0.16" @@ -42,9 +50,3 @@ TiledIteration = "0.5" TransformsBase = "1.6" Unitful = "1.17" julia = "1.9" - -[extensions] -MeshesMakieExt = "Makie" - -[weakdeps] -Makie = "ee78f7c6-11fb-53f2-987a-cfe4a2b5a57a" diff --git a/src/Meshes.jl b/src/Meshes.jl index db2612395..1f2f6efcb 100644 --- a/src/Meshes.jl +++ b/src/Meshes.jl @@ -4,6 +4,7 @@ module Meshes +using BinaryTrees using CoordRefSystems using StaticArrays using SparseArrays diff --git a/src/utils.jl b/src/utils.jl index adc3a611f..1a74be7cd 100644 --- a/src/utils.jl +++ b/src/utils.jl @@ -11,3 +11,4 @@ include("utils/cmp.jl") include("utils/units.jl") include("utils/crs.jl") include("utils/misc.jl") +include("utils/intersect.jl") diff --git a/src/utils/intersect.jl b/src/utils/intersect.jl new file mode 100644 index 000000000..4911028f5 --- /dev/null +++ b/src/utils/intersect.jl @@ -0,0 +1,236 @@ +# ------------------------------------------------------------------ +# Licensed under the MIT License. See LICENSE in the project root. +# ------------------------------------------------------------------ + +""" + bentleyottmann(segments; [digits]) + +Compute pairwise intersections between n `segments` +with `digits` precision in O(nβ‹…log(n)) time using +Bentley-Ottmann sweep line algorithm. + +By default, set `digits` based on the absolute +tolerance of the length type of the segments. + +## References + +* Bentley & Ottmann 1979. [Algorithms for reporting and counting + geometric intersections](https://ieeexplore.ieee.org/document/1675432) +""" +function bentleyottmann(segments; digits=_digits(segments)) + # orient segments + segs = map(segments) do s + a, b = extrema(s) + a > b ? (b, a) : (a, b) + end + + # retrieve types + P = eltype(first(segs)) + S = Tuple{P,P} + + # initialization + 𝒬 = BinaryTrees.AVLTree{P}() + β„› = BinaryTrees.AVLTree{S}() + ℬ = Dict{P,Vector{S}}() + β„° = Dict{P,Vector{S}}() + β„³ = Dict{P,Vector{S}}() + lookup = Dict{S,Int}() + for (i, (a, b)) in enumerate(segs) + BinaryTrees.insert!(𝒬, a) + BinaryTrees.insert!(𝒬, b) + haskey(ℬ, a) ? push!(ℬ[a], (a, b)) : (ℬ[a] = [(a, b)]) + haskey(β„°, b) ? push!(β„°[b], (a, b)) : (β„°[b] = [(a, b)]) + lookup[(a, b)] = i + end + + # sweep line algorithm + points = Vector{P}() + seginds = Vector{Vector{Int}}() + while !BinaryTrees.isempty(𝒬) + # current point (or event) + p = BinaryTrees.key(BinaryTrees.minnode(𝒬)) + + # delete point from event queue + BinaryTrees.delete!(𝒬, p) + # handle event, i.e. update 𝒬, β„› and β„³ + β„¬β‚š = get(ℬ, p, S[]) # segments with p at the begin + β„°β‚š = get(β„°, p, S[]) # segments with p at the end + β„³β‚š = get(β„³, p, S[]) # segments with p at the middle + # handle status line + _handlestatus!(β„›, β„¬β‚š, β„³β‚š, β„°β‚š) + + if length(β„³β‚š) > 0 + push!(points, p) + inds = [lookup[s] for s in β„³β‚š] + push!(seginds, inds) + end + + activesegs = β„¬β‚š βˆͺ β„³β‚š + + if isempty(activesegs) + for s in β„°β‚š + sβ‚—, sα΅£ = BinaryTrees.prevnext(β„›, s) + isnothing(sβ‚—) || isnothing(sα΅£) || Meshes._newevent!(𝒬, β„³, p, BinaryTrees.key(sβ‚—), BinaryTrees.key(sα΅£), digits) + end + else + _handlebottom!(activesegs[1], β„›, 𝒬, β„³, p, digits) + + _handletop!(activesegs[end], β„›, 𝒬, β„³, p, digits) + end + end + + points, seginds +end + +function _handlestatus!(β„›, β„¬β‚š, β„³β‚š, β„°β‚š) + for s in β„°β‚š + !isnothing(BinaryTrees.search(β„›, s)) || deleteat!(β„°β‚š, findfirst(isequal(s), β„°β‚š)) || BinaryTrees.delete!(β„›, s) + end + + for s in β„³β‚š + !isnothing(BinaryTrees.search(β„›, s)) || deleteat!(β„³β‚š, findfirst(isequal(s), β„³β‚š)) || BinaryTrees.delete!(β„›, s) + end + + for s in β„¬β‚š βˆͺ β„³β‚š + BinaryTrees.insert!(β„›, s) + end +end + +function _handlebottom!(s, β„›, 𝒬, β„³, p, digits) + sβ€² = BinaryTrees.key(BinaryTrees.search(β„›, s)) + + sβ‚—, _ = BinaryTrees.prevnext(β„›, sβ€²) + if !isnothing(sβ‚—) + _newevent!(𝒬, β„³, p, BinaryTrees.key(sβ‚—), sβ€², digits) + end +end + +function _handletop!(s, β„›, 𝒬, β„³, p, digits) + sβ€²β€² = BinaryTrees.search(β„›, s) |> BinaryTrees.key + + if !isnothing(sβ€²β€²) + _, sα΅£ = BinaryTrees.prevnext(β„›, sβ€²β€²) + if !isnothing(sα΅£) + _newevent!(𝒬, β„³, p, sβ€²β€², BinaryTrees.key(sα΅£), digits) + end + end +end + +# #sort by where segment intersects plane then check to β„›? #doesn't seem to work +# function _sort(segs, p) +# segs = copy(segs) +# if isempty(segs) +# return segs +# end +# _sort!(segs, p) +# end +# function _sort!(segs, p) +# T = lentype(eltype(first(segs))) +# ys = Vector{T}() +# pβ‚“, _ = coords(p) |> CoordRefSystems.values +# for s in segs +# a, b = s +# x₁, y₁ = coords(a) |> CoordRefSystems.values +# xβ‚‚, yβ‚‚ = coords(b) |> CoordRefSystems.values +# if x₁ == xβ‚‚ +# push!(ys, yβ‚‚) # Vertical segment, use yβ‚‚ directly to place at end +# else +# t = (pβ‚“ - x₁) / (xβ‚‚ - x₁) +# c = y₁ + t * (yβ‚‚ - y₁) +# push!(ys, c) +# end +# end +# segs[sortperm(ys)] +# end + +# _handlebeg!(β„¬β‚š, 𝒬, β„›, β„³, digits) +# _handleend!(β„°β‚š, 𝒬, β„›, β„³, digits) +# _handlemid!(β„³β‚š, 𝒬, β„›, β„³, digits) + +# report intersection point and segment indices +# inds = [lookup[s] for s in β„³β‚š] +# if !isempty(inds) +# if p ∈ keys(visited) +# seginds[visited[p]] = inds +# else +# push!(points, p) +# push!(seginds, inds) +# push!(visited, p => counter) +# counter += 1 +# end +# end +# function _handlebeg!(β„¬β‚š, 𝒬, β„›, β„³, digits) +# for s in β„¬β‚š +# BinaryTrees.insert!(β„›, s) +# prev, next = BinaryTrees.prevnext(β„›, s) +# isnothing(prev) || _newevent!(𝒬, β„³, BinaryTrees.key(prev), s, digits) +# isnothing(next) || _newevent!(𝒬, β„³, s, BinaryTrees.key(next), digits) +# isnothing(prev) || isnothing(next) || _rmevent!(𝒬, s, s, digits) +# end +# end + +# function _handleend!(β„°β‚š, 𝒬, β„›, β„³, digits) +# for s in β„°β‚š +# prev, next = BinaryTrees.prevnext(β„›, s) +# isnothing(prev) || isnothing(next) || _newevent!(𝒬, β„³, BinaryTrees.key(prev), BinaryTrees.key(next), digits) +# BinaryTrees.delete!(β„›, s) +# end +# end + +# function _handlemid!(β„³β‚š, 𝒬, β„›, β„³, digits) +# for s in β„³β‚š +# prev, next = BinaryTrees.prevnext(β„›, s) +# r = !isnothing(prev) ? BinaryTrees.key(prev) : nothing +# t = !isnothing(next) ? BinaryTrees.key(next) : nothing +# if !isnothing(r) +# _newevent!(𝒬, β„³, r, s, digits) +# if !isnothing(t) +# _newevent!(𝒬, β„³, r, t, digits) +# end +# end +# if !isnothing(t) +# _, next = BinaryTrees.prevnext(β„›, BinaryTrees.key(next)) +# u = !isnothing(next) ? BinaryTrees.key(next) : nothing +# if !isnothing(u) +# _newevent!(𝒬, β„³, t, u, digits) +# if !isnothing(r) +# _newevent!(𝒬, β„³, r, u, digits) +# end +# end +# end +# end +# end + +function _newevent!(𝒬, β„³, p, s₁, sβ‚‚, digits) + intersection(Segment(s₁), Segment(sβ‚‚)) do I + if type(I) == Crossing || type(I) == EdgeTouching + pβ€² = roundcoords(get(I); digits) + if pβ€² > p + !isnothing(BinaryTrees.search(𝒬, pβ€²)) || BinaryTrees.insert!(𝒬, pβ€²) + if haskey(β„³, pβ€²) + push!(β„³[pβ€²], s₁, sβ‚‚) + else + β„³[pβ€²] = [s₁, sβ‚‚] + end + end + end + nothing + end +end + +function _rmevent!(𝒬, s₁, sβ‚‚, digits) + intersection(Segment(s₁), Segment(sβ‚‚)) do I + if type(I) == Crossing || type(I) == EdgeTouching + p = roundcoords(get(I); digits) + BinaryTrees.delete!(𝒬, p) + end + nothing + end +end + +function _digits(segments) + s = first(segments) + β„’ = lentype(s) + Ο„ = ustrip(atol(β„’)) + round(Int, -log10(Ο„)) - 1 +end diff --git a/test/utils.jl b/test/utils.jl index 1db9fa283..9cb2d85e0 100644 --- a/test/utils.jl +++ b/test/utils.jl @@ -62,3 +62,63 @@ @test Meshes.roundcoords(pβ‚‚, digits=10) == p₁ @inferred Meshes.roundcoords(p₁, digits=10) end + +@testitem "bentleyottmann" setup = [Setup] begin + # basic check with a small number of segments + segs = + Segment.([ + (cart(0, 0), cart(1.1, 1.1)), + (cart(1, 0), cart(0, 1)), + (cart(0, 0), cart(0, 1)), + (cart(0, 0), cart(1, 0)), + (cart(0, 1), cart(1, 1)), + (cart(1, 0), cart(1, 1)) + ]) + points, seginds = Meshes.bentleyottmann(segs) + @test all(points .β‰ˆ [cart(0.5, 0.5), cart(1, 1)]) + @test length(points) == 2 + @test length(seginds) == 2 + @test seginds == [[1, 2], [1, 5, 6]] + @inferred Meshes.bentleyottmann(segs) + + segs = + Segment.([ + (cart(9, 13), cart(6, 9)), + (cart(2, 12), cart(9, 4.8)), + (cart(12, 11), cart(4, 7)), + (cart(2.5, 10), cart(12.5, 2)), + (cart(13, 6), cart(10, 4)), + (cart(10.5, 5.5), cart(9, 1)), + (cart(10, 4), cart(11, -1)), + (cart(10, 3), cart(10, 5)) + ]) + points, seginds = Meshes.bentleyottmann(segs) + @test length(points) == 4 + @test length(seginds) == 4 + @test Set(reduce(vcat, seginds)) == Set(2:8) + @test points[findfirst(p -> p β‰ˆ cart(10, 4), points)] β‰ˆ cart(10, 4) + @test Set(seginds[findfirst(p -> p β‰ˆ cart(10, 4), points)]) == Set([4, 5, 6, 7, 8]) + @test Set(seginds[findfirst(p -> p β‰ˆ cart(9, 4.8), points)]) == Set([4, 2]) + + # finds all intersections in a grid + horizontal = [Segment((1, i), (n, i)) for i in 1:n] + vertical = [Segment((i, 1), (i, n)) for i in 1:n] + segs = [horizontal; vertical] + points, seginds = Meshes.bentleyottmann(segs) + @test length(points) == 121 + @test length(seginds) == 121 + @test Set(length.(seginds)) == Set([2]) + + # result is invariant under rotations + segs = collect(segs) + for ΞΈ in T(Ο€ / 6):T(Ο€ / 6):T(2Ο€ - Ο€ / 6) + ΞΈpoints, ΞΈseginds = Meshes.bentleyottmann(segs |> Rotate(ΞΈ)) + @test length(ΞΈpoints) == 121 + @test length(ΞΈseginds) == 121 + @test Set(length.(ΞΈseginds)) == Set([2]) + end + + # inference test + segs = facets(cartgrid(10, 10)) + @inferred Meshes.bentleyottmann(segs) +end