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

updated LQ sideband selection logic #124

Merged
merged 2 commits into from
Feb 28, 2025
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
50 changes: 40 additions & 10 deletions ext/LegendSpecFitsRecipesBaseExt.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1390,9 +1390,9 @@
xlabel := "Drift Time"
ylabel := "LQ (A.U.)"
framestyle := :box
left_margin := -2Plots.mm
bottom_margin := -4Plots.mm
top_margin := -3Plots.mm
left_margin := (-2, :mm)
bottom_margin := (-4, :mm)
top_margin := (-3, :mm)

Check warning on line 1395 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1393-L1395

Added lines #L1393 - L1395 were not covered by tests
color := :viridis
formatter := :plain
thickness_scaling := 1.6
Expand Down Expand Up @@ -1476,15 +1476,15 @@

# recipe for the lq_cut report

@recipe function f(report::NamedTuple{(:cut, :fit_result, :temp_hists, :fit_report)}, lq_class::Vector{Float64}, e_cal, plot_type::Symbol)
@recipe function f(report::NamedTuple{(:cut, :fit_result, :temp_hists, :fit_report, :dep_σ, :edges)}, lq_class::Vector{Float64}, e_cal, plot_type::Symbol)

Check warning on line 1479 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1479

Added line #L1479 was not covered by tests

# Extract cutvalue from the report
cut_value = Measurements.value.(report.cut)

# Plot configuration for all types
left_margin := -2Plots.mm
bottom_margin := -4Plots.mm
top_margin := -3Plots.mm
left_margin := (-2, :mm)
bottom_margin := (-4, :mm)
top_margin := (-3, :mm)

Check warning on line 1487 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1485-L1487

Added lines #L1485 - L1487 were not covered by tests
thickness_scaling := 1.6
size := (1200, 900)
framestyle := :box
Expand Down Expand Up @@ -1662,10 +1662,40 @@
report.temp_hists.hist_sb2
end


elseif plot_type == :energy_spectrum

Check warning on line 1665 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1665

Added line #L1665 was not covered by tests
# Plot energy spectrum with DEP and sideband regions
left_margin := (-2, :mm)
bottom_margin := (-4, :mm)
top_margin := (-3, :mm)
thickness_scaling := 1.6
size := (1200, 900)
framestyle := :box
formatter := :plain
xlabel := "Energy"
ylabel := "Counts"

Check warning on line 1675 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1667-L1675

Added lines #L1667 - L1675 were not covered by tests

@series begin
seriestype := :stephist
label := "Energy Spectrum (σ: $(round(u"keV", report.dep_σ, sigdigits=3)))"
bins := 1000
e_cal[1500u"keV" .< e_cal .< 1660u"keV"]

Check warning on line 1681 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1677-L1681

Added lines #L1677 - L1681 were not covered by tests
end

@series begin
seriestype := :vline
label := "DEP Region"
fillcolor := :red
[report.edges.DEP_edge_left, report.edges.DEP_edge_right]

Check warning on line 1688 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1684-L1688

Added lines #L1684 - L1688 were not covered by tests
end

@series begin
seriestype := :vline
label := "Sideband Edges"
fillcolor := :green
legend := :topleft
[report.edges.sb1_edge, report.edges.sb2_edge]

Check warning on line 1696 in ext/LegendSpecFitsRecipesBaseExt.jl

View check run for this annotation

Codecov / codecov/patch

ext/LegendSpecFitsRecipesBaseExt.jl#L1691-L1696

Added lines #L1691 - L1696 were not covered by tests
end
end
end



end # module LegendSpecFitsRecipesBaseExt
30 changes: 25 additions & 5 deletions src/lqcut.jl
Original file line number Diff line number Diff line change
Expand Up @@ -159,11 +159,29 @@
dep_µ::Unitful.Energy, dep_σ::Unitful.Energy, e_cal::Vector{<:Unitful.Energy}, lq_classifier::Vector{<:AbstractFloat}; cut_sigma::Float64=3.0, dep_sideband_sigma::Float64=4.5, cut_truncation_sigma::Float64=3.5, uncertainty::Bool=true
)

# define sidebands
lq_dep = lq_classifier[dep_µ - dep_sideband_sigma * dep_σ .< e_cal .< dep_µ + dep_sideband_sigma * dep_σ]
lq_sb1 = lq_classifier[dep_µ - 2 * dep_sideband_sigma * dep_σ .< e_cal .< dep_µ - dep_sideband_sigma * dep_σ]
lq_sb2 = lq_classifier[dep_µ + dep_sideband_sigma * dep_σ .< e_cal .< dep_µ + 2 * dep_sideband_sigma * dep_σ]

# define sidebands; different for low and high energy resolution detectors to avoid sb reaching into 212-Bi FEP
DEP_edge_left = dep_µ - dep_sideband_sigma * dep_σ
DEP_edge_right = dep_µ + dep_sideband_sigma * dep_σ

if dep_σ < 2.0u"keV"
sb1_edge = dep_µ - 2 * dep_sideband_sigma * dep_σ
sb2_edge = dep_µ + 2 * dep_sideband_sigma * dep_σ

lq_dep = lq_classifier[DEP_edge_left .< e_cal .< DEP_edge_right]
lq_sb1 = lq_classifier[sb1_edge .< e_cal .< DEP_edge_left]
lq_sb2 = lq_classifier[DEP_edge_right .< e_cal .< sb2_edge]
else
sb1_edge = dep_µ - 2 * dep_sideband_sigma * dep_σ
sb2_edge = dep_µ - 3 * dep_sideband_sigma * dep_σ

Check warning on line 175 in src/lqcut.jl

View check run for this annotation

Codecov / codecov/patch

src/lqcut.jl#L174-L175

Added lines #L174 - L175 were not covered by tests

lq_dep = lq_classifier[DEP_edge_left .< e_cal .< DEP_edge_right]
lq_sb1 = lq_classifier[sb1_edge .< e_cal .< DEP_edge_left]
lq_sb2 = lq_classifier[sb2_edge .< e_cal .< sb1_edge]

Check warning on line 179 in src/lqcut.jl

View check run for this annotation

Codecov / codecov/patch

src/lqcut.jl#L177-L179

Added lines #L177 - L179 were not covered by tests
end

# save edges for crosschecks
edges_for_crosschecks = (;DEP_edge_left, DEP_edge_right, sb1_edge, sb2_edge)

# generate values for histogram edges
combined = filter(isfinite,[lq_dep; lq_sb1; lq_sb2])
ideal_bin_width = get_friedman_diaconis_bin_width(combined)
Expand Down Expand Up @@ -206,6 +224,8 @@
fit_result = fit_result,
temp_hists = temp_hists,
fit_report = fit_report,
dep_σ = dep_σ,
edges = edges_for_crosschecks,
)

return result, report
Expand Down
Loading