diff --git a/src/charts.jl b/src/charts.jl index b74ead2..beb761d 100644 --- a/src/charts.jl +++ b/src/charts.jl @@ -217,7 +217,7 @@ function _normals(tangents::SVector{2,SVector{2,T}}, ::Type{Val{0}}) where {T} t = tangents[1] s = tangents[2] - v = (t[1]*s[2] - t[2]*s[1])/2 + v = abs(t[1]*s[2] - t[2]*s[1])/2 # n[3] = tangents[1] × tangents[2] # l = norm(n) diff --git a/src/overlap.jl b/src/overlap.jl index bef6432..ef94ed3 100644 --- a/src/overlap.jl +++ b/src/overlap.jl @@ -106,3 +106,35 @@ function overlap(p::Simplex{3,3,0,4,T}, q::Simplex{3,3,0,4,T}) where T all(0+tol .<= u .<= 1-tol) && return true return false end + + + +function overlap(p::Simplex{2,2,0,3,T}, q::Simplex{2,2,0,3,T}) where T + + # tol = sqrt(eps(T)) + tol = 1e3 * eps(T) + + # Are the patches in the same plane? + u1 = q.tangents[1] + u2 = q.tangents[2] + v = p.vertices[1] - q.vertices[2] + + for i in 1:3 + a = p.vertices[mod1(i+1,3)] + b = p.vertices[mod1(i+2,3)] + c = p.vertices[i] + t = b - a + m = StaticArrays.SVector{2,T}(t[2],-t[1]) + + sp = zeros(T,3); sp[i] = dot(c-a, m) + sq = T[dot(q.vertices[j]-a, m) for j in 1:3] + + minp, maxp = extrema(sp) + minq, maxq = extrema(sq) + + maxq <= minp + tol && return false + maxp <= minq + tol && return false + end + + return true +end diff --git a/src/plotlyjs_glue.jl b/src/plotlyjs_glue.jl index 9fc87d2..d57f0e3 100644 --- a/src/plotlyjs_glue.jl +++ b/src/plotlyjs_glue.jl @@ -44,6 +44,44 @@ function __init__() return s end + @eval function patch(Γ::CompScienceMeshes.AbstractMesh{2}, fcr=nothing; + caxis=nothing, showscale=true, color="red", kwargs...) + + v = vertexarray(Γ) + c = cellarray(Γ) + + x = v[:,1]; y = v[:,2]; z = zeros(length(y)) + i = c[:,1].-1; j = c[:,2].-1; k = c[:,3].-1 + + + if fcr != nothing + m, M = extrema(fcr) + if caxis != nothing + m, M = caxis + end + + s = PlotlyBase.mesh3d(; + x=x, y=y, z=z, + i=i, j=j, k=k, + intensitymode="cell", + intensity=fcr, + colorscale="Viridis", + showscale=showscale, + cmin=m, + cmax=M, + kwargs... + ) + else + s = PlotlyBase.mesh3d(; + x=x, y=y, z=z, + i=i, j=j, k=k, + color=color, + kwargs... + ) + end + return s + end + @eval function patch(a::Vector{<:Simplex}; kwargs...) vertices = reduce(vcat, [v.vertices for v in a]) faces = collect(SVector(3*(i-1)+1, 3*(i-1)+2, 3*(i-1)+3) for i in 1:length(a)) diff --git a/test/test_overlap.jl b/test/test_overlap.jl index 437fca0..b1c34d2 100644 --- a/test/test_overlap.jl +++ b/test/test_overlap.jl @@ -2,59 +2,92 @@ using Test using CompScienceMeshes using StaticArrays -p = simplex( - point(0,0,0), - point(1,0,0), - point(0,1,0)) - -q1 = simplex( - point(0.6, 0.6, 0.0), - point(1.6, 0.6, 0.0), - point(0.6, 1.6, 0.0), -) - -q2 = simplex( - point(0.4, 0.4, 0.0), - point(1.4, 0.4, 0.0), - point(0.4, 1.4, 0.0), -) - -@test overlap(p, q1) == false -@test overlap(p, q2) == true - -q1 = simplex( - point(0.949726718617726,-0.278864925637586,0.142314838273285), - point(0.989821441880933,0.0,-0.142314838273285), - point(0.989821441880933,0.0,0.142314838273285), -) - -q2 = simplex( - point(0.949726718617726,-0.278864925637586,-0.142314838273285), - point(0.989821441880933,0.0,-0.142314838273285), - point(0.949726718617726,-0.278864925637586,0.142314838273285), -) - -@test overlap(q1, q2) == false - -## make sure the submesh function work for 1D meshes - - -l1 = meshsegment(1.0,1/2) -vt = skeleton(l1,0) -bd = boundary(l1) - -overlaps = overlap_gpredicate(bd) -# pred1 = c -> overlaps(simplex(vertices(vt,c))) -pred1 = c -> overlaps(chart(vt,c)) -@test pred1(CompScienceMeshes.SimplexGraph(1)) -@test !pred1(CompScienceMeshes.SimplexGraph(2)) -@test pred1(CompScienceMeshes.SimplexGraph(3)) - - -# test a case where the segments are: -# not of unit length -# colinear and opposite -# meet in a common point -ch1 = simplex(point(1/3,0,0), point(1/3,1/3,0)) -ch2 = simplex(point(1/3,1/3,0), point(1/3,2/3,0)) -@test !overlap(ch1, ch2) +let # udim 3 + p = simplex( + point(0,0,0), + point(1,0,0), + point(0,1,0)) + + q1 = simplex( + point(0.6, 0.6, 0.0), + point(1.6, 0.6, 0.0), + point(0.6, 1.6, 0.0), + ) + + q2 = simplex( + point(0.4, 0.4, 0.0), + point(1.4, 0.4, 0.0), + point(0.4, 1.4, 0.0), + ) + + @test overlap(p, q1) == false + @test overlap(p, q2) == true + + q1 = simplex( + point(0.949726718617726,-0.278864925637586,0.142314838273285), + point(0.989821441880933,0.0,-0.142314838273285), + point(0.989821441880933,0.0,0.142314838273285), + ) + + q2 = simplex( + point(0.949726718617726,-0.278864925637586,-0.142314838273285), + point(0.989821441880933,0.0,-0.142314838273285), + point(0.949726718617726,-0.278864925637586,0.142314838273285), + ) + + @test overlap(q1, q2) == false + + ## make sure the submesh function work for 1D meshes + + + l1 = meshsegment(1.0,1/2) + vt = skeleton(l1,0) + bd = boundary(l1) + + overlaps = overlap_gpredicate(bd) + # pred1 = c -> overlaps(simplex(vertices(vt,c))) + pred1 = c -> overlaps(chart(vt,c)) + @test pred1(CompScienceMeshes.SimplexGraph(1)) + @test !pred1(CompScienceMeshes.SimplexGraph(2)) + @test pred1(CompScienceMeshes.SimplexGraph(3)) + + + # test a case where the segments are: + # not of unit length + # colinear and opposite + # meet in a common point + ch1 = simplex(point(1/3,0,0), point(1/3,1/3,0)) + ch2 = simplex(point(1/3,1/3,0), point(1/3,2/3,0)) + @test !overlap(ch1, ch2) +end + +let # udim 2 + p = simplex( + point(0,0), + point(1,0), + point(0,1)) + + q1 = simplex( + point(0.6, 0.6), + point(1.6, 0.6), + point(0.6, 1.6), + ) + + q2 = simplex( + point(0.4, 0.4), + point(1.4, 0.4), + point(0.4, 1.4), + ) + + @test overlap(p, q1) == false + @test overlap(p, q2) == true + + + # test a case where the segments are: + # not of unit length + # colinear and opposite + # meet in a common point + ch1 = simplex(point(1/3,0), point(1/3,1/3)) + ch2 = simplex(point(1/3,1/3), point(1/3,2/3)) + @test !overlap(ch1, ch2) +end \ No newline at end of file diff --git a/test/test_patches.jl b/test/test_patches.jl index 403c156..1dbbebc 100644 --- a/test/test_patches.jl +++ b/test/test_patches.jl @@ -1,6 +1,7 @@ using Test using CompScienceMeshes +# universe dimension 3 for T in [Float32, Float64] local mesh = meshrectangle(T(1.0), T(1.0), T(1.0)) # local faces = skeleton(mesh, 2) @@ -19,3 +20,23 @@ for T in [Float32, Float64] point(T, 0.0, 0.0, -1.0)] @test p.volume == T(0.5) end + + +# universe dimension 2 +for T in [Float32, Float64] + local mesh = meshrectangle(T(1.0), T(1.0), T(1.0),2) + # local faces = skeleton(mesh, 2) + # local verts = cellvertices(mesh, 1) + p = chart(mesh, cells(mesh)[1]) + #p = simplex(verts) + + @test p.vertices == [ + point(T, 0.0, 0.0), + point(T, 0.0, 1.0), + point(T, 1.0, 0.0)] + @test p.tangents == [ + point(T, -1.0, 0.0), + point(T, -1.0, 1.0)] + @test p.normals == [] + @test p.volume == T(0.5) +end