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

Mesh Plotting #218

Merged
merged 15 commits into from
Oct 14, 2021
10 changes: 2 additions & 8 deletions src/ConstructiveSolidGeometry/SurfacePrimitives/ConeMantle.jl
Original file line number Diff line number Diff line change
Expand Up @@ -60,15 +60,9 @@ function normal(cm::ConeMantle{T,Tuple{T,T},<:Any,:outwards}, pt::CartesianPoint
CartesianPoint(CylindricalPoint{T}( one(T), cyl.φ, -Δr / Δz)), cm))
end

function _get_n_points_in_arc(cm::ConeMantle, n::Int64)::Int64
φMin, φMax = get_φ_limits(cm)
f = (φMax - φMin)/(2π)
Int(ceil(n*f))
end

function vertices(cm::ConeMantle{T}, n::Int64)::Vector{CartesianPoint{T}} where {T}
φMin, φMax = get_φ_limits(cm)
n = _get_n_points_in_arc(cm, n)
n = _get_n_points_in_arc_φ(cm, n)
φ = range(φMin, φMax, length = n+1)
rbot = radius_at_z(cm,-cm.hZ)
rtop = radius_at_z(cm,cm.hZ)
Expand All @@ -78,7 +72,7 @@ function vertices(cm::ConeMantle{T}, n::Int64)::Vector{CartesianPoint{T}} where
end

function connections(cm::ConeMantle, n::Int64)::Vector{Vector{Int64}}
n = _get_n_points_in_arc(cm, n)
n = _get_n_points_in_arc_φ(cm, n)
[[i,i+1,i+n+2,i+n+1] for i in 1:n]
end

Expand Down
12 changes: 0 additions & 12 deletions src/ConstructiveSolidGeometry/SurfacePrimitives/EllipsoidMantle.jl
Original file line number Diff line number Diff line change
Expand Up @@ -76,18 +76,6 @@ function normal(em::EllipsoidMantle{T,T,TP,TT,:outwards}, pt::CartesianPoint{T})
end
normal(em::EllipsoidMantle{T,TR,TP,TT,:inwards}, pt::CartesianPoint{T}) where {T,TR,TP,TT} = -normal(flip(em), pt)

function _get_n_points_in_arc_φ(em::EllipsoidMantle, n::Int64)::Int64
φMin, φMax = get_φ_limits(em)
f = (φMax - φMin)/(2π)
Int(ceil(n*f))
end

function _get_n_points_in_arc_θ(em::EllipsoidMantle, n::Int64)::Int64
θMin, θMax = get_θ_limits(em)
f = (θMax - θMin)/(π)
Int(ceil(n*f))
end

function vertices(em::EllipsoidMantle{T}, n::Int64)::Vector{CartesianPoint{T}} where {T}
rx, ry, rz = get_radii(em)
φMin, φMax = get_φ_limits(em)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -36,34 +36,28 @@ Plane(es::EllipticalSurface{T}) where {T} = Plane{T}(es.origin, es.rotation * Ca

normal(es::EllipticalSurface{T}, ::CartesianPoint{T} = zero(CartesianPoint{T})) where {T} = es.rotation * CartesianVector{T}(zero(T), zero(T), one(T))

function _get_n_points_in_arc(es::EllipticalSurface, n::Int64)::Int64
φMin, φMax = get_φ_limits(es)
f = (φMax - φMin)/(2π)
Int(ceil(n*f))
end

function vertices(es::EllipticalSurface{T, T}, n::Int64)::Vector{CartesianPoint{T}} where {T}
φMin, φMax = get_φ_limits(es)
n = _get_n_points_in_arc(es, n)
n = _get_n_points_in_arc_φ(es, n)
φ = range(φMin, φMax, length = n+1)
append!([es.origin],[_transform_into_global_coordinate_system(CartesianPoint{T}(es.r*cos(φ), es.r*sin(φ), 0), es) for φ in φ])
end

function vertices(es::EllipticalSurface{T, Tuple{T,T}}, n::Int64)::Vector{CartesianPoint{T}} where {T}
rMin, rMax = es.r
φMin, φMax = get_φ_limits(es)
n = _get_n_points_in_arc(es, n)
n = _get_n_points_in_arc_φ(es, n)
φ = range(φMin, φMax, length = n+1)
[_transform_into_global_coordinate_system(CartesianPoint{T}(r*cos(φ), r*sin(φ), 0), es) for r in (rMin,rMax) for φ in φ]
end

function connections(es::EllipticalSurface{T, T}, n::Int64)::Vector{Vector{Int64}} where {T}
n = _get_n_points_in_arc(es, n)
n = _get_n_points_in_arc_φ(es, n)
[[1,i,i+1] for i in 2:n+1]
end

function connections(es::EllipticalSurface{T, Tuple{T,T}}, n::Int64)::Vector{Vector{Int64}} where {T}
n = _get_n_points_in_arc(es, n)
n = _get_n_points_in_arc_φ(es, n)
[[i,i+1,i+n+2,i+n+1] for i in 1:n]
end

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -11,3 +11,15 @@ include("TorusMantle.jl")
# Normal vectors of planar surface primitives do not depend on the point,
# but the normal method is generally called with a point.
normal(p::AbstractPlanarSurfacePrimitive, ::AbstractCoordinatePoint) = normal(p)

function _get_n_points_in_arc_φ(p::AbstractSurfacePrimitive, n::Int64)::Int64
φMin, φMax = get_φ_limits(p)
f = (φMax - φMin)/(2π)
Int(ceil(n*f))
end

function _get_n_points_in_arc_θ(p::AbstractSurfacePrimitive, n::Int64)::Int64
θMin, θMax = get_θ_limits(p)
f = (θMax - θMin)/(2π)
Int(ceil(n*f))
end
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I am getting an error that these functions are not defined when trying to plot a detector.
I guess they need to go before all the includes in this file.

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Seems to be on your side. I was just able to plot a detector locally.
Also, the docs build went through where a detector is plotted.
I used:

using Plots, SolidStateDetectors
T = Float32 
sim = Simulation{T}(SSD_examples[:InvertedCoax])
plot(sim.detector)

12 changes: 0 additions & 12 deletions src/ConstructiveSolidGeometry/SurfacePrimitives/TorusMantle.jl
Original file line number Diff line number Diff line change
Expand Up @@ -46,18 +46,6 @@ function normal(tm::TorusMantle{T,TP,TT,:outwards}, pt::CartesianPoint{T}) where
end
normal(tm::TorusMantle{T,TP,TT,:inwards}, pt::CartesianPoint{T}) where {T,TP,TT} = -normal(flip(tm), pt)

function _get_n_points_in_arc_φ(tm::TorusMantle, n::Int64)::Int64
φMin, φMax = get_φ_limits(tm)
f = (φMax - φMin)/(2π)
Int(ceil(n*f))
end

function _get_n_points_in_arc_θ(tm::TorusMantle, n::Int64)::Int64
θMin, θMax = get_θ_limits(tm)
f = (θMax - θMin)/(2π)
Int(ceil(n*f))
end

function vertices(tm::TorusMantle{T}, n::Int64)::Vector{CartesianPoint{T}} where {T}
φMin, φMax = get_φ_limits(tm)
θMin, θMax = get_θ_limits(tm)
Expand Down
8 changes: 8 additions & 0 deletions src/ConstructiveSolidGeometry/plotting/CSG/CSG.jl
Original file line number Diff line number Diff line change
Expand Up @@ -10,12 +10,20 @@ end
ps = primitives(csg)
@series begin
label --> "CSG"
if !isClosedPrimitive(ps[1])
fillcolor := :white
fillalpha --> 0.2
end
linestyle := isClosedPrimitive(ps[1]) ? :solid : :dot
ps[1]
end
for i in 2:length(ps)
@series begin
label := nothing
if !isClosedPrimitive(ps[i])
fillcolor := :white
fillalpha --> 0.2
end
linestyle := isClosedPrimitive(ps[i]) ? :solid : :dot
ps[i]
end
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,6 @@
if haskey(plotattributes, :seriestype) && plotattributes[:seriestype] == :mesh3d
@series begin
label --> "Cone Mantle"
println(n)
mesh(cm, n = n)
end
else
Expand Down
2 changes: 2 additions & 0 deletions src/PlotRecipes/SolidStateDetector.jl
Original file line number Diff line number Diff line change
Expand Up @@ -12,6 +12,8 @@ end

@recipe function f(contact::Contact)
linecolor --> contact.id
fillcolor --> contact.id
fillalpha --> 0.2
l = contact.name != "" ? "$(contact.name) (id: $(contact.id))" : "Contact - id: $(contact.id)"
label --> l
contact.geometry
Expand Down