diff --git a/src/pseudo_prior.jl b/src/pseudo_prior.jl index 4ef0936f..cdcbc151 100644 --- a/src/pseudo_prior.jl +++ b/src/pseudo_prior.jl @@ -3,29 +3,36 @@ function get_standard_pseudo_prior(h::Histogram, ps::NamedTuple{(:peak_pos, :pea μ = ifelse(fixed_position, ConstValueDist(ps.peak_pos), Normal(ps.peak_pos, 0.2*ps.peak_sigma)), σ = weibull_from_mx(ps.peak_sigma, 1.5*ps.peak_sigma), n = weibull_from_mx(ps.peak_counts, 1.5*ps.peak_counts), - step_amplitude = weibull_from_mx(ps.mean_background_step, ps.mean_background_step + 5*ps.mean_background_std), skew_fraction = ifelse(low_e_tail, truncated(weibull_from_mx(0.002, 0.008), 0.0, 0.5), ConstValueDist(0.0)), skew_width = ifelse(low_e_tail, weibull_from_mx(ps.peak_sigma/ps.peak_pos, 2*ps.peak_sigma/ps.peak_pos), ConstValueDist(1.0)), background = weibull_from_mx(ps.mean_background, ps.mean_background + 5*ps.mean_background_std), ) - if fit_func == :f_fit + pprior_step = NamedTupleDist(step_amplitude = weibull_from_mx(ps.mean_background_step, ps.mean_background_step + 5*ps.mean_background_std),) + if fit_func == :gamma_def + return merge(pprior_base, pprior_step) + elseif fit_func == :gamma_bckFlat return pprior_base + elseif fit_func == :gamma_tails || fit_func == :gamma_tails_bckFlat + pprior_tails = merge(pprior_base, NamedTupleDist( + skew_fraction_highE = ifelse(low_e_tail, truncated(weibull_from_mx(0.01, 0.05), 0.0, 0.1), ConstValueDist(0.0)), + skew_width_highE = ifelse(low_e_tail, weibull_from_mx(0.001, 1e-2), ConstValueDist(1.0)),) + ) + if fit_func == :gamma_tails + return merge(pprior_tails, pprior_step) + else + return pprior_tails + end else window_left = ps.peak_pos - minimum(h.edges[1]) window_right = maximum(h.edges[1]) - ps.peak_pos - return merge(pprior_base, - if fit_func == :f_fit_bckSlope + return merge(pprior_base, pprior_step, + if fit_func == :gamma_bckSlope NamedTupleDist( background_slope = ifelse(ps.mean_background < 5, ConstValueDist(0), truncated(Normal(0, 0.1*ps.mean_background_std / (window_left + window_right)), - ps.mean_background / window_right, 0)), ) - elseif fit_func ==:f_fit_bckExp + elseif fit_func ==:gamma_bckExp NamedTupleDist(background_exp = weibull_from_mx(3e-2, 5e-2)) - elseif fit_func == :f_fit_tails - NamedTupleDist( - skew_fraction_highE = ifelse(low_e_tail, truncated(weibull_from_mx(0.01, 0.05), 0.0, 0.1), ConstValueDist(0.0)), - skew_width_highE = ifelse(low_e_tail, weibull_from_mx(0.001, 1e-2), ConstValueDist(1.0)), - ) end ) end diff --git a/src/specfit.jl b/src/specfit.jl index ea0a10f7..cfcb99b4 100644 --- a/src/specfit.jl +++ b/src/specfit.jl @@ -22,7 +22,8 @@ end export fit_peaks function fit_peaks_th228(peakhists::Array, peakstats::StructArray, th228_lines::Vector{T},; e_unit::Union{Nothing, Unitful.EnergyUnits}=nothing, uncertainty::Bool=true, low_e_tail::Bool=true, iterative_fit::Bool=false, - fit_func::Symbol= :f_fit, pseudo_prior::NamedTupleDist=NamedTupleDist(empty = true), m_cal_simple::MaybeWithEnergyUnits = 1.0) where T<:Any + + fit_func::Vector{Symbol}= fill(:gamma_def, length(th228_lines)), pseudo_prior::NamedTupleDist=NamedTupleDist(empty = true), m_cal_simple::MaybeWithEnergyUnits = 1.0) where T<:Any e_unit = ifelse(isnothing(e_unit), NoUnits, e_unit) @@ -30,22 +31,23 @@ function fit_peaks_th228(peakhists::Array, peakstats::StructArray, th228_lines:: v_result = Vector{NamedTuple}(undef, length(th228_lines)) v_report = Vector{NamedTuple}(undef, length(th228_lines)) - # iterate throuh all peaks Threads.@threads for i in eachindex(th228_lines) # get histogram and peakstats peak = th228_lines[i] h = peakhists[i] ps = peakstats[i] + f = fit_func[i] + # fit peak - result_peak, report_peak = fit_single_peak_th228(h, ps; uncertainty=uncertainty, low_e_tail=low_e_tail, fit_func = fit_func, pseudo_prior = pseudo_prior) + result_peak, report_peak = fit_single_peak_th228(h, ps; uncertainty=uncertainty, low_e_tail=low_e_tail, fit_func = f, pseudo_prior = pseudo_prior) # check covariance matrix for being semi positive definite (no negative uncertainties) if uncertainty if iterative_fit && !isposdef(result_peak.covmat) @warn "Covariance matrix not positive definite for peak $peak - repeat fit without low energy tail" pval_save = result_peak.pval - result_peak, report_peak = fit_single_peak_th228(h, ps, ; uncertainty=uncertainty, low_e_tail=false, fit_func = fit_func, pseudo_prior = pseudo_prior) + result_peak, report_peak = fit_single_peak_th228(h, ps, ; uncertainty=uncertainty, low_e_tail=false, fit_func = f, pseudo_prior = pseudo_prior) @info "New covariance matrix is positive definite: $(isposdef(result_peak.covmat))" @info "p-val with low-energy tail p=$(round(pval_save,digits=5)) , without low-energy tail: p=$(round((result_peak.pval),digits=5))" end @@ -76,7 +78,7 @@ Also, FWHM is calculated from the fitted peakshape with MC error propagation. Th """ function fit_single_peak_th228(h::Histogram, ps::NamedTuple{(:peak_pos, :peak_fwhm, :peak_sigma, :peak_counts, :bin_width, :mean_background, :mean_background_step, :mean_background_std), NTuple{8, T}}; uncertainty::Bool=true, low_e_tail::Bool=true, fixed_position::Bool=false, pseudo_prior::NamedTupleDist=NamedTupleDist(empty = true), - fit_func::Symbol=:f_fit, background_center::Union{Real,Nothing} = ps.peak_pos, m_cal_simple::Real = 1.0) where T<:Real + fit_func::Symbol=:gamma_def, background_center::Union{Real,Nothing} = ps.peak_pos, m_cal_simple::Real = 1.0) where T<:Real # create standard pseudo priors pseudo_prior = get_pseudo_prior(h, ps, fit_func; pseudo_prior = pseudo_prior, fixed_position = fixed_position, low_e_tail = low_e_tail) @@ -142,7 +144,8 @@ function fit_single_peak_th228(h::Histogram, ps::NamedTuple{(:peak_pos, :peak_fw @debug "FWHM: $(fwhm) ± $(fwhm_err)" result = merge(NamedTuple{keys(v_ml)}([measurement(v_ml[k], v_ml_err[k]) for k in keys(v_ml)]...), - (fwhm = measurement(fwhm, fwhm_err), + (fwhm = measurement(fwhm, fwhm_err), + fit_func = fit_func, gof = (pvalue = pval, chi2 = chi2, dof = dof, @@ -170,7 +173,7 @@ function fit_single_peak_th228(h::Histogram, ps::NamedTuple{(:peak_pos, :peak_fw @debug "FWHM: $(fwhm)" result = merge(NamedTuple{keys(v_ml)}([measurement(v_ml[k], NaN) for k in keys(v_ml)]...), - (fwhm = measurement(fwhm, NaN), ), (gof = (converged = converged,),)) + (fwhm = measurement(fwhm, NaN), fit_func = fit_func, ), (gof = (converged = converged,),)) report = ( v = v_ml, h = h, @@ -208,7 +211,7 @@ Get the FWHM of a peak from the fit parameters. """ function estimate_fwhm(v::NamedTuple) # get FWHM - f_sigWithTail = Base.Fix2(get_th228_fit_functions().f_sigWithTail,v) + f_sigWithTail = Base.Fix2(get_th228_fit_functions().gamma_sigWithTail,v) try if v.skew_fraction <= 0.5 half_max_sig = maximum(f_sigWithTail.(v.μ - v.σ:0.001:v.μ + v.σ))/2 @@ -268,7 +271,7 @@ Also, FWHM is calculated from the fitted peakshape with MC error propagation. Th function fit_subpeaks_th228( h_survived::Histogram, h_cut::Histogram, h_result; uncertainty::Bool=false, low_e_tail::Bool=true, fix_σ::Bool = true, fix_skew_fraction::Bool = true, fix_skew_width::Bool = true, - pseudo_prior::NamedTupleDist=NamedTupleDist(empty = true), fit_func::Symbol=:f_fit, background_center::Real = h_result.μ + pseudo_prior::NamedTupleDist=NamedTupleDist(empty = true), fit_func::Symbol= :gamma_def, background_center::Real = h_result.μ ) # create standard pseudo priors @@ -437,7 +440,7 @@ function fit_subpeaks_th228( @debug "σ cut : $(v_ml.σ_cut)" result = merge(NamedTuple{keys(v_ml)}([measurement(v_ml[k], NaN) for k in keys(v_ml)]...), - (fwhm = measurement(fwhm, NaN), gof = (converged = converged, survived = NamedTuple(), cut = NamedTuple()))) + (fwhm = measurement(fwhm, NaN), fit_func = fit_func, gof = (converged = converged, survived = NamedTuple(), cut = NamedTuple()))) end report = ( diff --git a/src/specfit_functions.jl b/src/specfit_functions.jl index ab403906..0bf0f7ff 100644 --- a/src/specfit_functions.jl +++ b/src/specfit_functions.jl @@ -1,27 +1,29 @@ # This file is a part of LegendSpecFits.jl, licensed under the MIT License (MIT). - -# helper functions for fitting peakshapes, legacy -th228_fit_functions = ( - f_fit = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background), - f_fit_bckSlope = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_slope = v.background_slope), - f_fit_bckExp = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_exp = v.background_exp), - f_fit_tails = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width , v.background; skew_fraction_highE = v.skew_fraction_highE, skew_width_highE = v.skew_width_highE), - f_sigWithTail = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction) + lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), -) - +""" + get_th228_fit_functions(; background_center::Union{Real,Nothing} = nothing) +This function defines the gamma peakshape fit functions used in the calibration specfits. +* gamma_def: "default" gamma peakshape with gaussian signal, low-energy tail, and background (flat + step) +* gamma_tails: default gamma peakshape + high-energy tail +* gamma_bckSlope: default gamma peakshape + linear background slope +* gamma_bckExp: default gamma peakshape + exponential background +* gamma_bckFlat: default gamma peakshape - step background (only flat component!) +* gamma_tails_bckFlat: default gamma peakshape + high-energy tail - step background (only flat component!) +""" function get_th228_fit_functions(; background_center::Union{Real,Nothing} = nothing) merge( - (f_fit = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background), - f_fit_tails = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width , v.background; skew_fraction_highE = v.skew_fraction_highE, skew_width_highE = v.skew_width_highE), - f_sigWithTail = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction) + lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), + (gamma_def = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background), + gamma_tails = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width , v.background; skew_fraction_highE = v.skew_fraction_highE, skew_width_highE = v.skew_width_highE), + gamma_sigWithTail = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction) + lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), + gamma_bckFlat = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, 0.0, v.skew_fraction, v.skew_width, v.background), + gamma_tails_bckFlat = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, 0.0, v.skew_fraction, v.skew_width , v.background; skew_fraction_highE = v.skew_fraction_highE, skew_width_highE = v.skew_width_highE), ), if isnothing(background_center) - (f_fit_bckSlope = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_slope = v.background_slope, background_center = v.μ), - f_fit_bckExp = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_exp = v.background_exp, background_center = v.μ), + (gamma_bckSlope = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_slope = v.background_slope, background_center = v.μ), + gamma_bckExp = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_exp = v.background_exp, background_center = v.μ), ) else - (f_fit_bckSlope = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_slope = v.background_slope, background_center = background_center), - f_fit_bckExp = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_exp = v.background_exp, background_center = background_center), + (gamma_bckSlope = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_slope = v.background_slope, background_center = background_center), + gamma_bckExp = (x, v) -> gamma_peakshape(x, v.μ, v.σ, v.n, v.step_amplitude, v.skew_fraction, v.skew_width, v.background; background_exp = v.background_exp, background_center = background_center), ) end ) @@ -33,23 +35,32 @@ This function defines the components (signal, low/high-energy tail, backgrounds) These component functions are used in the fit-report and in plot receipes """ function peakshape_components(fit_func::Symbol; background_center::Union{Real,Nothing} = nothing) - if fit_func == :f_fit + if fit_func == :gamma_def funcs = (f_sig = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction), f_lowEtail = (x, v) -> lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), f_bck = (x, v) -> background_peakshape(x, v.μ, v.σ, v.step_amplitude, v.background)) - elseif fit_func == :f_fit_bckSlope + elseif fit_func == :gamma_bckSlope funcs = (f_sig = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction), f_lowEtail = (x, v) -> lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), f_bck = (x, v) -> background_peakshape(x, v.μ, v.σ, v.step_amplitude, v.background; background_slope = v.background_slope, background_center = background_center)) - elseif fit_func == :f_fit_bckExp + elseif fit_func == :gamma_bckExp funcs = (f_sig = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction), f_lowEtail = (x, v) -> lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), f_bck = (x, v) -> background_peakshape(x, v.μ, v.σ, v.step_amplitude, v.background; background_exp = v.background_exp, background_center = background_center)) - elseif fit_func == :f_fit_tails + elseif fit_func == :gamma_bckFlat + funcs = (f_sig = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction), + f_lowEtail = (x, v) -> lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), + f_bck = (x, v) -> background_peakshape(x, v.μ, v.σ, 0.0, v.background)) + elseif fit_func == :gamma_tails funcs = (f_sig = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, (v.skew_fraction + v.skew_fraction_highE)), f_lowEtail = (x, v) -> lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), f_highEtail = (x, v) -> highEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction_highE, v.skew_width_highE), f_bck = (x, v) -> background_peakshape(x, v.μ, v.σ, v.step_amplitude, v.background)) + elseif fit_func == :gamma_tails_bckFlat + funcs = (f_sig = (x, v) -> signal_peakshape(x, v.μ, v.σ, v.n, (v.skew_fraction + v.skew_fraction_highE)), + f_lowEtail = (x, v) -> lowEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction, v.skew_width), + f_highEtail = (x, v) -> highEtail_peakshape(x, v.μ, v.σ, v.n, v.skew_fraction_highE, v.skew_width_highE), + f_bck = (x, v) -> background_peakshape(x, v.μ, v.σ, 0.0, v.background)) end labels = (f_sig = "Signal", f_lowEtail = "Low-energy tail", f_bck = "Background", f_highEtail = "High-energy tail") colors = (f_sig = :orangered1, f_lowEtail = :orange, f_bck = :dodgerblue2, f_highEtail = :forestgreen) diff --git a/test/test_utils.jl b/test/test_utils.jl index c1119416..cd918bb2 100644 --- a/test/test_utils.jl +++ b/test/test_utils.jl @@ -10,7 +10,7 @@ Sample Legend200 calibration data based on "Inverse Transform Sampling" method * find the value x such that F(x) = u by solving for x . done by interpolation of the inverse cdf * repeat for many u : energy samples """ -function generate_mc_spectrum(n_tot::Int=200000,; f_fit::Base.Callable=LegendSpecFits.get_th228_fit_functions().f_fit) +function generate_mc_spectrum(n_tot::Int=200000,; f_fit::Base.Callable=LegendSpecFits.get_th228_fit_functions().gamma_def) th228_lines = [583.191, 727.330, 860.564, 1592.53, 1620.50, 2103.53, 2614.51] @@ -53,7 +53,7 @@ function generate_mc_spectrum(n_tot::Int=200000,; f_fit::Base.Callable=LegendSpe for i=1:length(th228_lines) bandwidth = maximum(model_cdf_all[i])-minimum(model_cdf_all[i]) rand_i = minimum(model_cdf_all[i]).+bandwidth.*rand(n_i[i]); # make sure sample is within model range - interp_cdf_inv = linear_interpolation(model_cdf_all[i],bin_centers_all[i]) # inverse cdf + interp_cdf_inv = linear_interpolation(model_cdf_all[i], bin_centers_all[i]) # inverse cdf energy_mc_all[i] = interp_cdf_inv.(rand_i) end